1 2 #include <private/tsimpl.h> /*I "petscts.h" I*/ 3 4 /* Logging support */ 5 PetscClassId TS_CLASSID; 6 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "TSSetTypeFromOptions" 10 /* 11 TSSetTypeFromOptions - Sets the type of ts from user options. 12 13 Collective on TS 14 15 Input Parameter: 16 . ts - The ts 17 18 Level: intermediate 19 20 .keywords: TS, set, options, database, type 21 .seealso: TSSetFromOptions(), TSSetType() 22 */ 23 static PetscErrorCode TSSetTypeFromOptions(TS ts) 24 { 25 PetscBool opt; 26 const char *defaultType; 27 char typeName[256]; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 if (((PetscObject)ts)->type_name) { 32 defaultType = ((PetscObject)ts)->type_name; 33 } else { 34 defaultType = TSEULER; 35 } 36 37 if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 38 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 39 if (opt) { 40 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 41 } else { 42 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 43 } 44 PetscFunctionReturn(0); 45 } 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "TSSetFromOptions" 49 /*@ 50 TSSetFromOptions - Sets various TS parameters from user options. 51 52 Collective on TS 53 54 Input Parameter: 55 . ts - the TS context obtained from TSCreate() 56 57 Options Database Keys: 58 + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 59 . -ts_max_steps maxsteps - maximum number of time-steps to take 60 . -ts_max_time time - maximum time to compute to 61 . -ts_dt dt - initial time step 62 . -ts_monitor - print information at each timestep 63 - -ts_monitor_draw - plot information at each timestep 64 65 Level: beginner 66 67 .keywords: TS, timestep, set, options, database 68 69 .seealso: TSGetType() 70 @*/ 71 PetscErrorCode TSSetFromOptions(TS ts) 72 { 73 PetscReal dt; 74 PetscBool opt,flg; 75 PetscErrorCode ierr; 76 PetscViewer monviewer; 77 char monfilename[PETSC_MAX_PATH_LEN]; 78 SNES snes; 79 80 PetscFunctionBegin; 81 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 82 ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr); 83 84 /* Handle generic TS options */ 85 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 88 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 89 if (opt) { 90 ts->initial_time_step = ts->time_step = dt; 91 } 92 93 /* Monitor options */ 94 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 95 if (flg) { 96 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 97 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 98 } 99 opt = PETSC_FALSE; 100 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 101 if (opt) { 102 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 103 } 104 opt = PETSC_FALSE; 105 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 106 if (opt) { 107 ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 108 } 109 110 /* Handle TS type options */ 111 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 112 113 /* Handle specific TS options */ 114 if (ts->ops->setfromoptions) { 115 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 116 } 117 118 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 119 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 120 ierr = PetscOptionsEnd();CHKERRQ(ierr); 121 122 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 123 /* Handle subobject options */ 124 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 125 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "TSViewFromOptions" 131 /*@ 132 TSViewFromOptions - This function visualizes the ts based upon user options. 133 134 Collective on TS 135 136 Input Parameter: 137 . ts - The ts 138 139 Level: intermediate 140 141 .keywords: TS, view, options, database 142 .seealso: TSSetFromOptions(), TSView() 143 @*/ 144 PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 145 { 146 PetscViewer viewer; 147 PetscDraw draw; 148 PetscBool opt = PETSC_FALSE; 149 char fileName[PETSC_MAX_PATH_LEN]; 150 PetscErrorCode ierr; 151 152 PetscFunctionBegin; 153 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 154 if (opt && !PetscPreLoadingOn) { 155 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 156 ierr = TSView(ts, viewer);CHKERRQ(ierr); 157 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 158 } 159 opt = PETSC_FALSE; 160 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 161 if (opt) { 162 ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 163 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 164 if (title) { 165 ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 166 } else { 167 ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 168 ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 169 } 170 ierr = TSView(ts, viewer);CHKERRQ(ierr); 171 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 172 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 173 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 174 } 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "TSComputeRHSJacobian" 180 /*@ 181 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 182 set with TSSetRHSJacobian(). 183 184 Collective on TS and Vec 185 186 Input Parameters: 187 + ts - the TS context 188 . t - current timestep 189 - x - input vector 190 191 Output Parameters: 192 + A - Jacobian matrix 193 . B - optional preconditioning matrix 194 - flag - flag indicating matrix structure 195 196 Notes: 197 Most users should not need to explicitly call this routine, as it 198 is used internally within the nonlinear solvers. 199 200 See KSPSetOperators() for important information about setting the 201 flag parameter. 202 203 TSComputeJacobian() is valid only for TS_NONLINEAR 204 205 Level: developer 206 207 .keywords: SNES, compute, Jacobian, matrix 208 209 .seealso: TSSetRHSJacobian(), KSPSetOperators() 210 @*/ 211 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 212 { 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 217 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 218 PetscCheckSameComm(ts,1,X,3); 219 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only"); 220 if (ts->ops->rhsjacobian) { 221 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 222 *flg = DIFFERENT_NONZERO_PATTERN; 223 PetscStackPush("TS user Jacobian function"); 224 ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 225 PetscStackPop; 226 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 227 /* make sure user returned a correct Jacobian and preconditioner */ 228 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 229 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 230 } else { 231 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 232 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 233 if (*A != *B) { 234 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 236 } 237 *flg = ts->rhsmatstructure; 238 } 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "TSComputeRHSFunction" 244 /*@ 245 TSComputeRHSFunction - Evaluates the right-hand-side function. 246 247 Collective on TS and Vec 248 249 Input Parameters: 250 + ts - the TS context 251 . t - current time 252 - x - state vector 253 254 Output Parameter: 255 . y - right hand side 256 257 Note: 258 Most users should not need to explicitly call this routine, as it 259 is used internally within the nonlinear solvers. 260 261 If the user did not provide a function but merely a matrix, 262 this routine applies the matrix. 263 264 Level: developer 265 266 .keywords: TS, compute 267 268 .seealso: TSSetRHSFunction(), TSComputeIFunction() 269 @*/ 270 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 271 { 272 PetscErrorCode ierr; 273 274 PetscFunctionBegin; 275 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 276 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 277 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 278 279 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 280 if (ts->ops->rhsfunction) { 281 PetscStackPush("TS user right-hand-side function"); 282 ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 283 PetscStackPop; 284 } else if (ts->ops->rhsmatrix) { 285 Mat A,B; 286 ierr = TSGetMatrices(ts,&A,&B,PETSC_NULL, PETSC_NULL,PETSC_NULL,PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 287 ts->rhsmatstructure = DIFFERENT_NONZERO_PATTERN; 288 PetscStackPush("TS user right-hand-side matrix function"); 289 ierr = (*ts->ops->rhsmatrix)(ts,t,&A,&B,&ts->rhsmatstructure,ts->jacP);CHKERRQ(ierr); 290 PetscStackPop; 291 /* call TSSetMatrices() in case the user changed the pointers */ 292 ierr = TSSetMatrices(ts,A,B,PETSC_NULL, PETSC_NULL,PETSC_NULL,PETSC_NULL, PETSC_NULL,ts->jacP);CHKERRQ(ierr); 293 ierr = MatMult(A,x,y);CHKERRQ(ierr); 294 } else SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"No RHS provided, must call TSSetRHSFunction() or TSSetMatrices()"); 295 296 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 297 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "TSComputeIFunction" 303 /*@ 304 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 305 306 Collective on TS and Vec 307 308 Input Parameters: 309 + ts - the TS context 310 . t - current time 311 . X - state vector 312 - Xdot - time derivative of state vector 313 314 Output Parameter: 315 . Y - right hand side 316 317 Note: 318 Most users should not need to explicitly call this routine, as it 319 is used internally within the nonlinear solvers. 320 321 If the user did did not write their equations in implicit form, this 322 function recasts them in implicit form. 323 324 Level: developer 325 326 .keywords: TS, compute 327 328 .seealso: TSSetIFunction(), TSComputeRHSFunction() 329 @*/ 330 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y) 331 { 332 PetscErrorCode ierr; 333 334 PetscFunctionBegin; 335 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 336 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 337 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 338 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 339 340 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 341 if (ts->ops->ifunction) { 342 PetscStackPush("TS user implicit function"); 343 ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 344 PetscStackPop; 345 } else { 346 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 347 /* Convert to implicit form: F(X,Xdot) = Alhs * Xdot - Frhs(X) */ 348 if (ts->ops->lhsmatrix) { 349 ts->lhsmatstructure = DIFFERENT_NONZERO_PATTERN; 350 PetscStackPush("TS user left-hand-side matrix function"); 351 ierr = (*ts->ops->lhsmatrix)(ts,t,&ts->Alhs,&ts->Blhs,&ts->lhsmatstructure,ts->jacP);CHKERRQ(ierr); 352 PetscStackPop; 353 ierr = VecScale(Y,-1.);CHKERRQ(ierr); 354 ierr = MatMultAdd(ts->Alhs,Xdot,Y,Y);CHKERRQ(ierr); 355 } else { 356 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 357 } 358 } 359 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "TSComputeIJacobian" 365 /*@ 366 TSComputeIJacobian - Evaluates the Jacobian of the DAE 367 368 Collective on TS and Vec 369 370 Input 371 Input Parameters: 372 + ts - the TS context 373 . t - current timestep 374 . X - state vector 375 . Xdot - time derivative of state vector 376 - shift - shift to apply, see note below 377 378 Output Parameters: 379 + A - Jacobian matrix 380 . B - optional preconditioning matrix 381 - flag - flag indicating matrix structure 382 383 Notes: 384 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 385 386 dF/dX + shift*dF/dXdot 387 388 Most users should not need to explicitly call this routine, as it 389 is used internally within the nonlinear solvers. 390 391 TSComputeIJacobian() is valid only for TS_NONLINEAR 392 393 Level: developer 394 395 .keywords: TS, compute, Jacobian, matrix 396 397 .seealso: TSSetIJacobian() 398 @*/ 399 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg) 400 { 401 PetscErrorCode ierr; 402 403 PetscFunctionBegin; 404 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 405 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 406 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 407 PetscValidPointer(A,6); 408 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 409 PetscValidPointer(B,7); 410 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 411 PetscValidPointer(flg,8); 412 413 *flg = SAME_NONZERO_PATTERN; /* In case it we're solving a linear problem in which case it wouldn't get initialized below. */ 414 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 415 if (ts->ops->ijacobian) { 416 PetscStackPush("TS user implicit Jacobian"); 417 ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 418 PetscStackPop; 419 } else { 420 if (ts->ops->rhsjacobian) { 421 PetscStackPush("TS user right-hand-side Jacobian"); 422 ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 423 PetscStackPop; 424 } else { 425 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 426 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 427 if (*A != *B) { 428 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 429 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 430 } 431 } 432 433 /* Convert to implicit form */ 434 /* inefficient because these operations will normally traverse all matrix elements twice */ 435 ierr = MatScale(*A,-1);CHKERRQ(ierr); 436 if (ts->Alhs) { 437 ierr = MatAXPY(*A,shift,ts->Alhs,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 438 } else { 439 ierr = MatShift(*A,shift);CHKERRQ(ierr); 440 } 441 if (*A != *B) { 442 ierr = MatScale(*B,-1);CHKERRQ(ierr); 443 ierr = MatShift(*B,shift);CHKERRQ(ierr); 444 } 445 } 446 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNCT__ 451 #define __FUNCT__ "TSSetRHSFunction" 452 /*@C 453 TSSetRHSFunction - Sets the routine for evaluating the function, 454 F(t,u), where U_t = F(t,u). 455 456 Logically Collective on TS 457 458 Input Parameters: 459 + ts - the TS context obtained from TSCreate() 460 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 461 . f - routine for evaluating the right-hand-side function 462 - ctx - [optional] user-defined context for private data for the 463 function evaluation routine (may be PETSC_NULL) 464 465 Calling sequence of func: 466 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 467 468 + t - current timestep 469 . u - input vector 470 . F - function vector 471 - ctx - [optional] user-defined function context 472 473 Important: 474 The user MUST call either this routine or TSSetMatrices(). 475 476 Level: beginner 477 478 .keywords: TS, timestep, set, right-hand-side, function 479 480 .seealso: TSSetMatrices() 481 @*/ 482 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 483 { 484 PetscErrorCode ierr; 485 SNES snes; 486 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 489 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 490 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 491 ts->ops->rhsfunction = f; 492 ts->funP = ctx; 493 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 494 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "TSSetMatrices" 500 /*@C 501 TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs, 502 where Alhs(t) U_t = Arhs(t) U. 503 504 Logically Collective on TS 505 506 Input Parameters: 507 + ts - the TS context obtained from TSCreate() 508 . Arhs - matrix 509 . Brhs - preconditioning matrix 510 . frhs - the matrix evaluation routine for Arhs and Brhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 511 if Arhs is not a function of t. 512 . Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix. 513 . Blhs - preconditioning matrix or PETSC_NULL if Blhs is an identity matrix 514 . flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 515 if Alhs is not a function of t. 516 . flag - flag indicating information about the matrix structure of Arhs and Alhs. 517 The available options are 518 SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs 519 DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs 520 - ctx - [optional] user-defined context for private data for the 521 matrix evaluation routine (may be PETSC_NULL) 522 523 Calling sequence of func: 524 $ func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 525 526 + t - current timestep 527 . A - matrix A, where U_t = A(t) U 528 . B - preconditioner matrix, usually the same as A 529 . flag - flag indicating information about the preconditioner matrix 530 structure (same as flag in KSPSetOperators()) 531 - ctx - [optional] user-defined context for matrix evaluation routine 532 533 Notes: 534 The routine func() takes Mat* as the matrix arguments rather than Mat. 535 This allows the matrix evaluation routine to replace Arhs or Alhs with a 536 completely new new matrix structure (not just different matrix elements) 537 when appropriate, for instance, if the nonzero structure is changing 538 throughout the global iterations. 539 540 Important: 541 The user MUST call either this routine or TSSetRHSFunction(). 542 543 Level: beginner 544 545 .keywords: TS, timestep, set, matrix 546 547 .seealso: TSSetRHSFunction() 548 @*/ 549 PetscErrorCode TSSetMatrices(TS ts,Mat Arhs,Mat Brhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,Mat Blhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx) 550 { 551 PetscErrorCode ierr; 552 SNES snes; 553 554 PetscFunctionBegin; 555 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 556 if (frhs) ts->ops->rhsmatrix = frhs; 557 if (flhs) ts->ops->lhsmatrix = flhs; 558 if (ctx) ts->jacP = ctx; 559 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 560 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 561 ierr = SNESSetJacobian(snes,Arhs,Brhs,SNESTSFormJacobian,ts);CHKERRQ(ierr); 562 if (Alhs) { 563 PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4); 564 PetscCheckSameComm(ts,1,Alhs,4); 565 ierr = PetscObjectReference((PetscObject)Alhs);CHKERRQ(ierr); 566 ierr = MatDestroy(&ts->Alhs);CHKERRQ(ierr); 567 ts->Alhs = Alhs; 568 } 569 if (Blhs) { 570 PetscValidHeaderSpecific(Blhs,MAT_CLASSID,4); 571 PetscCheckSameComm(ts,1,Blhs,4); 572 ierr = PetscObjectReference((PetscObject)Blhs);CHKERRQ(ierr); 573 ierr = MatDestroy(&ts->Blhs);CHKERRQ(ierr); 574 ts->Blhs = Blhs; 575 } 576 ts->rhsmatstructure = flag; 577 ts->lhsmatstructure = flag; 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "TSGetMatrices" 583 /*@C 584 TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 585 where Alhs(t) U_t = Arhs(t) U. 586 587 Not Collective, but parallel objects are returned if TS is parallel 588 589 Input Parameter: 590 . ts - The TS context obtained from TSCreate() 591 592 Output Parameters: 593 + Arhs - The right-hand side matrix 594 . Alhs - The left-hand side matrix 595 - ctx - User-defined context for matrix evaluation routine 596 597 Notes: You can pass in PETSC_NULL for any return argument you do not need. 598 599 Level: intermediate 600 601 .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 602 603 .keywords: TS, timestep, get, matrix 604 605 @*/ 606 PetscErrorCode TSGetMatrices(TS ts,Mat *Arhs,Mat *Brhs,TSMatrix *frhs,Mat *Alhs,Mat *Blhs,TSMatrix *flhs,void **ctx) 607 { 608 PetscErrorCode ierr; 609 SNES snes; 610 611 PetscFunctionBegin; 612 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 613 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 614 ierr = SNESGetJacobian(snes,Arhs,Brhs,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 615 if (frhs) *frhs = ts->ops->rhsmatrix; 616 if (Alhs) *Alhs = ts->Alhs; 617 if (Blhs) *Blhs = ts->Blhs; 618 if (flhs) *flhs = ts->ops->lhsmatrix; 619 if (ctx) *ctx = ts->jacP; 620 PetscFunctionReturn(0); 621 } 622 623 #undef __FUNCT__ 624 #define __FUNCT__ "TSSetRHSJacobian" 625 /*@C 626 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 627 where U_t = F(U,t), as well as the location to store the matrix. 628 Use TSSetMatrices() for linear problems. 629 630 Logically Collective on TS 631 632 Input Parameters: 633 + ts - the TS context obtained from TSCreate() 634 . A - Jacobian matrix 635 . B - preconditioner matrix (usually same as A) 636 . f - the Jacobian evaluation routine 637 - ctx - [optional] user-defined context for private data for the 638 Jacobian evaluation routine (may be PETSC_NULL) 639 640 Calling sequence of func: 641 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 642 643 + t - current timestep 644 . u - input vector 645 . A - matrix A, where U_t = A(t)u 646 . B - preconditioner matrix, usually the same as A 647 . flag - flag indicating information about the preconditioner matrix 648 structure (same as flag in KSPSetOperators()) 649 - ctx - [optional] user-defined context for matrix evaluation routine 650 651 Notes: 652 See KSPSetOperators() for important information about setting the flag 653 output parameter in the routine func(). Be sure to read this information! 654 655 The routine func() takes Mat * as the matrix arguments rather than Mat. 656 This allows the matrix evaluation routine to replace A and/or B with a 657 completely new matrix structure (not just different matrix elements) 658 when appropriate, for instance, if the nonzero structure is changing 659 throughout the global iterations. 660 661 Level: beginner 662 663 .keywords: TS, timestep, set, right-hand-side, Jacobian 664 665 .seealso: TSDefaultComputeJacobianColor(), 666 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 667 668 @*/ 669 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 670 { 671 PetscErrorCode ierr; 672 SNES snes; 673 674 PetscFunctionBegin; 675 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 676 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 677 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 678 if (A) PetscCheckSameComm(ts,1,A,2); 679 if (B) PetscCheckSameComm(ts,1,B,3); 680 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 681 682 if (f) ts->ops->rhsjacobian = f; 683 if (ctx) ts->jacP = ctx; 684 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 685 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 686 PetscFunctionReturn(0); 687 } 688 689 690 #undef __FUNCT__ 691 #define __FUNCT__ "TSSetIFunction" 692 /*@C 693 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 694 695 Logically Collective on TS 696 697 Input Parameters: 698 + ts - the TS context obtained from TSCreate() 699 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 700 . f - the function evaluation routine 701 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 702 703 Calling sequence of f: 704 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 705 706 + t - time at step/stage being solved 707 . u - state vector 708 . u_t - time derivative of state vector 709 . F - function vector 710 - ctx - [optional] user-defined context for matrix evaluation routine 711 712 Important: 713 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 714 715 Level: beginner 716 717 .keywords: TS, timestep, set, DAE, Jacobian 718 719 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 720 @*/ 721 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 722 { 723 PetscErrorCode ierr; 724 SNES snes; 725 726 PetscFunctionBegin; 727 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 728 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 729 if (f) ts->ops->ifunction = f; 730 if (ctx) ts->funP = ctx; 731 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 732 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "TSGetIFunction" 738 /*@C 739 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 740 741 Not Collective 742 743 Input Parameter: 744 . ts - the TS context 745 746 Output Parameter: 747 + r - vector to hold residual (or PETSC_NULL) 748 . func - the function to compute residual (or PETSC_NULL) 749 - ctx - the function context (or PETSC_NULL) 750 751 Level: advanced 752 753 .keywords: TS, nonlinear, get, function 754 755 .seealso: TSSetIFunction(), SNESGetFunction() 756 @*/ 757 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 758 { 759 PetscErrorCode ierr; 760 SNES snes; 761 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 764 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 765 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 766 if (func) *func = ts->ops->ifunction; 767 if (ctx) *ctx = ts->funP; 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "TSGetRHSFunction" 773 /*@C 774 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 775 776 Not Collective 777 778 Input Parameter: 779 . ts - the TS context 780 781 Output Parameter: 782 + r - vector to hold computed right hand side (or PETSC_NULL) 783 . func - the function to compute right hand side (or PETSC_NULL) 784 - ctx - the function context (or PETSC_NULL) 785 786 Level: advanced 787 788 .keywords: TS, nonlinear, get, function 789 790 .seealso: TSSetRhsfunction(), SNESGetFunction() 791 @*/ 792 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 793 { 794 PetscErrorCode ierr; 795 SNES snes; 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 799 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 800 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 801 if (func) *func = ts->ops->rhsfunction; 802 if (ctx) *ctx = ts->funP; 803 PetscFunctionReturn(0); 804 } 805 806 #undef __FUNCT__ 807 #define __FUNCT__ "TSSetIJacobian" 808 /*@C 809 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 810 you provided with TSSetIFunction(). 811 812 Logically Collective on TS 813 814 Input Parameters: 815 + ts - the TS context obtained from TSCreate() 816 . A - Jacobian matrix 817 . B - preconditioning matrix for A (may be same as A) 818 . f - the Jacobian evaluation routine 819 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 820 821 Calling sequence of f: 822 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 823 824 + t - time at step/stage being solved 825 . U - state vector 826 . U_t - time derivative of state vector 827 . a - shift 828 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 829 . B - preconditioning matrix for A, may be same as A 830 . flag - flag indicating information about the preconditioner matrix 831 structure (same as flag in KSPSetOperators()) 832 - ctx - [optional] user-defined context for matrix evaluation routine 833 834 Notes: 835 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 836 837 The matrix dF/dU + a*dF/dU_t you provide turns out to be 838 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 839 The time integrator internally approximates U_t by W+a*U where the positive "shift" 840 a and vector W depend on the integration method, step size, and past states. For example with 841 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 842 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 843 844 Level: beginner 845 846 .keywords: TS, timestep, DAE, Jacobian 847 848 .seealso: TSSetIFunction(), TSSetRHSJacobian() 849 850 @*/ 851 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 852 { 853 PetscErrorCode ierr; 854 SNES snes; 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 858 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 859 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 860 if (A) PetscCheckSameComm(ts,1,A,2); 861 if (B) PetscCheckSameComm(ts,1,B,3); 862 if (f) ts->ops->ijacobian = f; 863 if (ctx) ts->jacP = ctx; 864 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 865 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 866 PetscFunctionReturn(0); 867 } 868 869 #undef __FUNCT__ 870 #define __FUNCT__ "TSView" 871 /*@C 872 TSView - Prints the TS data structure. 873 874 Collective on TS 875 876 Input Parameters: 877 + ts - the TS context obtained from TSCreate() 878 - viewer - visualization context 879 880 Options Database Key: 881 . -ts_view - calls TSView() at end of TSStep() 882 883 Notes: 884 The available visualization contexts include 885 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 886 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 887 output where only the first processor opens 888 the file. All other processors send their 889 data to the first processor to print. 890 891 The user can open an alternative visualization context with 892 PetscViewerASCIIOpen() - output to a specified file. 893 894 Level: beginner 895 896 .keywords: TS, timestep, view 897 898 .seealso: PetscViewerASCIIOpen() 899 @*/ 900 PetscErrorCode TSView(TS ts,PetscViewer viewer) 901 { 902 PetscErrorCode ierr; 903 const TSType type; 904 PetscBool iascii,isstring,isundials; 905 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 908 if (!viewer) { 909 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 910 } 911 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 912 PetscCheckSameComm(ts,1,viewer,2); 913 914 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 915 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 916 if (iascii) { 917 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 918 if (ts->ops->view) { 919 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 920 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 921 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 922 } 923 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 924 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 925 if (ts->problem_type == TS_NONLINEAR) { 926 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 927 } 928 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 929 } else if (isstring) { 930 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 931 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 932 } 933 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 934 ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 935 if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 936 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 937 PetscFunctionReturn(0); 938 } 939 940 941 #undef __FUNCT__ 942 #define __FUNCT__ "TSSetApplicationContext" 943 /*@ 944 TSSetApplicationContext - Sets an optional user-defined context for 945 the timesteppers. 946 947 Logically Collective on TS 948 949 Input Parameters: 950 + ts - the TS context obtained from TSCreate() 951 - usrP - optional user context 952 953 Level: intermediate 954 955 .keywords: TS, timestep, set, application, context 956 957 .seealso: TSGetApplicationContext() 958 @*/ 959 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 960 { 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 963 ts->user = usrP; 964 PetscFunctionReturn(0); 965 } 966 967 #undef __FUNCT__ 968 #define __FUNCT__ "TSGetApplicationContext" 969 /*@ 970 TSGetApplicationContext - Gets the user-defined context for the 971 timestepper. 972 973 Not Collective 974 975 Input Parameter: 976 . ts - the TS context obtained from TSCreate() 977 978 Output Parameter: 979 . usrP - user context 980 981 Level: intermediate 982 983 .keywords: TS, timestep, get, application, context 984 985 .seealso: TSSetApplicationContext() 986 @*/ 987 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 988 { 989 PetscFunctionBegin; 990 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 991 *(void**)usrP = ts->user; 992 PetscFunctionReturn(0); 993 } 994 995 #undef __FUNCT__ 996 #define __FUNCT__ "TSGetTimeStepNumber" 997 /*@ 998 TSGetTimeStepNumber - Gets the current number of timesteps. 999 1000 Not Collective 1001 1002 Input Parameter: 1003 . ts - the TS context obtained from TSCreate() 1004 1005 Output Parameter: 1006 . iter - number steps so far 1007 1008 Level: intermediate 1009 1010 .keywords: TS, timestep, get, iteration, number 1011 @*/ 1012 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 1013 { 1014 PetscFunctionBegin; 1015 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1016 PetscValidIntPointer(iter,2); 1017 *iter = ts->steps; 1018 PetscFunctionReturn(0); 1019 } 1020 1021 #undef __FUNCT__ 1022 #define __FUNCT__ "TSSetInitialTimeStep" 1023 /*@ 1024 TSSetInitialTimeStep - Sets the initial timestep to be used, 1025 as well as the initial time. 1026 1027 Logically Collective on TS 1028 1029 Input Parameters: 1030 + ts - the TS context obtained from TSCreate() 1031 . initial_time - the initial time 1032 - time_step - the size of the timestep 1033 1034 Level: intermediate 1035 1036 .seealso: TSSetTimeStep(), TSGetTimeStep() 1037 1038 .keywords: TS, set, initial, timestep 1039 @*/ 1040 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 1041 { 1042 PetscFunctionBegin; 1043 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1044 ts->time_step = time_step; 1045 ts->initial_time_step = time_step; 1046 ts->ptime = initial_time; 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "TSSetTimeStep" 1052 /*@ 1053 TSSetTimeStep - Allows one to reset the timestep at any time, 1054 useful for simple pseudo-timestepping codes. 1055 1056 Logically Collective on TS 1057 1058 Input Parameters: 1059 + ts - the TS context obtained from TSCreate() 1060 - time_step - the size of the timestep 1061 1062 Level: intermediate 1063 1064 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1065 1066 .keywords: TS, set, timestep 1067 @*/ 1068 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1069 { 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1072 PetscValidLogicalCollectiveReal(ts,time_step,2); 1073 ts->time_step = time_step; 1074 PetscFunctionReturn(0); 1075 } 1076 1077 #undef __FUNCT__ 1078 #define __FUNCT__ "TSGetTimeStep" 1079 /*@ 1080 TSGetTimeStep - Gets the current timestep size. 1081 1082 Not Collective 1083 1084 Input Parameter: 1085 . ts - the TS context obtained from TSCreate() 1086 1087 Output Parameter: 1088 . dt - the current timestep size 1089 1090 Level: intermediate 1091 1092 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1093 1094 .keywords: TS, get, timestep 1095 @*/ 1096 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1097 { 1098 PetscFunctionBegin; 1099 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1100 PetscValidDoublePointer(dt,2); 1101 *dt = ts->time_step; 1102 PetscFunctionReturn(0); 1103 } 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "TSGetSolution" 1107 /*@ 1108 TSGetSolution - Returns the solution at the present timestep. It 1109 is valid to call this routine inside the function that you are evaluating 1110 in order to move to the new timestep. This vector not changed until 1111 the solution at the next timestep has been calculated. 1112 1113 Not Collective, but Vec returned is parallel if TS is parallel 1114 1115 Input Parameter: 1116 . ts - the TS context obtained from TSCreate() 1117 1118 Output Parameter: 1119 . v - the vector containing the solution 1120 1121 Level: intermediate 1122 1123 .seealso: TSGetTimeStep() 1124 1125 .keywords: TS, timestep, get, solution 1126 @*/ 1127 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1128 { 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1131 PetscValidPointer(v,2); 1132 *v = ts->vec_sol; 1133 PetscFunctionReturn(0); 1134 } 1135 1136 /* ----- Routines to initialize and destroy a timestepper ---- */ 1137 #undef __FUNCT__ 1138 #define __FUNCT__ "TSSetProblemType" 1139 /*@ 1140 TSSetProblemType - Sets the type of problem to be solved. 1141 1142 Not collective 1143 1144 Input Parameters: 1145 + ts - The TS 1146 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1147 .vb 1148 U_t = A U 1149 U_t = A(t) U 1150 U_t = F(t,U) 1151 .ve 1152 1153 Level: beginner 1154 1155 .keywords: TS, problem type 1156 .seealso: TSSetUp(), TSProblemType, TS 1157 @*/ 1158 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1159 { 1160 PetscFunctionBegin; 1161 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1162 ts->problem_type = type; 1163 PetscFunctionReturn(0); 1164 } 1165 1166 #undef __FUNCT__ 1167 #define __FUNCT__ "TSGetProblemType" 1168 /*@C 1169 TSGetProblemType - Gets the type of problem to be solved. 1170 1171 Not collective 1172 1173 Input Parameter: 1174 . ts - The TS 1175 1176 Output Parameter: 1177 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1178 .vb 1179 M U_t = A U 1180 M(t) U_t = A(t) U 1181 U_t = F(t,U) 1182 .ve 1183 1184 Level: beginner 1185 1186 .keywords: TS, problem type 1187 .seealso: TSSetUp(), TSProblemType, TS 1188 @*/ 1189 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1190 { 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1193 PetscValidIntPointer(type,2); 1194 *type = ts->problem_type; 1195 PetscFunctionReturn(0); 1196 } 1197 1198 #undef __FUNCT__ 1199 #define __FUNCT__ "TSSetUp" 1200 /*@ 1201 TSSetUp - Sets up the internal data structures for the later use 1202 of a timestepper. 1203 1204 Collective on TS 1205 1206 Input Parameter: 1207 . ts - the TS context obtained from TSCreate() 1208 1209 Notes: 1210 For basic use of the TS solvers the user need not explicitly call 1211 TSSetUp(), since these actions will automatically occur during 1212 the call to TSStep(). However, if one wishes to control this 1213 phase separately, TSSetUp() should be called after TSCreate() 1214 and optional routines of the form TSSetXXX(), but before TSStep(). 1215 1216 Level: advanced 1217 1218 .keywords: TS, timestep, setup 1219 1220 .seealso: TSCreate(), TSStep(), TSDestroy() 1221 @*/ 1222 PetscErrorCode TSSetUp(TS ts) 1223 { 1224 PetscErrorCode ierr; 1225 1226 PetscFunctionBegin; 1227 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1228 if (ts->setupcalled) PetscFunctionReturn(0); 1229 1230 if (!((PetscObject)ts)->type_name) { 1231 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1232 } 1233 1234 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1235 1236 if (ts->ops->setup) { 1237 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1238 } 1239 1240 ts->setupcalled = PETSC_TRUE; 1241 PetscFunctionReturn(0); 1242 } 1243 1244 #undef __FUNCT__ 1245 #define __FUNCT__ "TSReset" 1246 /*@ 1247 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1248 1249 Collective on TS 1250 1251 Input Parameter: 1252 . ts - the TS context obtained from TSCreate() 1253 1254 Level: beginner 1255 1256 .keywords: TS, timestep, reset 1257 1258 .seealso: TSCreate(), TSSetup(), TSDestroy() 1259 @*/ 1260 PetscErrorCode TSReset(TS ts) 1261 { 1262 PetscErrorCode ierr; 1263 1264 PetscFunctionBegin; 1265 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1266 if (ts->ops->reset) { 1267 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1268 } 1269 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1270 ierr = MatDestroy(&ts->Alhs);CHKERRQ(ierr); 1271 ierr = MatDestroy(&ts->Blhs);CHKERRQ(ierr); 1272 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1273 if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);} 1274 ts->setupcalled = PETSC_FALSE; 1275 PetscFunctionReturn(0); 1276 } 1277 1278 #undef __FUNCT__ 1279 #define __FUNCT__ "TSDestroy" 1280 /*@ 1281 TSDestroy - Destroys the timestepper context that was created 1282 with TSCreate(). 1283 1284 Collective on TS 1285 1286 Input Parameter: 1287 . ts - the TS context obtained from TSCreate() 1288 1289 Level: beginner 1290 1291 .keywords: TS, timestepper, destroy 1292 1293 .seealso: TSCreate(), TSSetUp(), TSSolve() 1294 @*/ 1295 PetscErrorCode TSDestroy(TS *ts) 1296 { 1297 PetscErrorCode ierr; 1298 1299 PetscFunctionBegin; 1300 if (!*ts) PetscFunctionReturn(0); 1301 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1302 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1303 1304 ierr = TSReset((*ts));CHKERRQ(ierr); 1305 1306 /* if memory was published with AMS then destroy it */ 1307 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1308 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1309 1310 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1311 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1312 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1313 1314 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1315 PetscFunctionReturn(0); 1316 } 1317 1318 #undef __FUNCT__ 1319 #define __FUNCT__ "TSGetSNES" 1320 /*@ 1321 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1322 a TS (timestepper) context. Valid only for nonlinear problems. 1323 1324 Not Collective, but SNES is parallel if TS is parallel 1325 1326 Input Parameter: 1327 . ts - the TS context obtained from TSCreate() 1328 1329 Output Parameter: 1330 . snes - the nonlinear solver context 1331 1332 Notes: 1333 The user can then directly manipulate the SNES context to set various 1334 options, etc. Likewise, the user can then extract and manipulate the 1335 KSP, KSP, and PC contexts as well. 1336 1337 TSGetSNES() does not work for integrators that do not use SNES; in 1338 this case TSGetSNES() returns PETSC_NULL in snes. 1339 1340 Level: beginner 1341 1342 .keywords: timestep, get, SNES 1343 @*/ 1344 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1345 { 1346 PetscErrorCode ierr; 1347 1348 PetscFunctionBegin; 1349 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1350 PetscValidPointer(snes,2); 1351 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1352 if (!ts->snes) { 1353 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1354 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1355 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1356 } 1357 *snes = ts->snes; 1358 PetscFunctionReturn(0); 1359 } 1360 1361 #undef __FUNCT__ 1362 #define __FUNCT__ "TSGetKSP" 1363 /*@ 1364 TSGetKSP - Returns the KSP (linear solver) associated with 1365 a TS (timestepper) context. 1366 1367 Not Collective, but KSP is parallel if TS is parallel 1368 1369 Input Parameter: 1370 . ts - the TS context obtained from TSCreate() 1371 1372 Output Parameter: 1373 . ksp - the nonlinear solver context 1374 1375 Notes: 1376 The user can then directly manipulate the KSP context to set various 1377 options, etc. Likewise, the user can then extract and manipulate the 1378 KSP and PC contexts as well. 1379 1380 TSGetKSP() does not work for integrators that do not use KSP; 1381 in this case TSGetKSP() returns PETSC_NULL in ksp. 1382 1383 Level: beginner 1384 1385 .keywords: timestep, get, KSP 1386 @*/ 1387 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1388 { 1389 PetscErrorCode ierr; 1390 SNES snes; 1391 1392 PetscFunctionBegin; 1393 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1394 PetscValidPointer(ksp,2); 1395 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1396 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1397 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1398 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1399 PetscFunctionReturn(0); 1400 } 1401 1402 /* ----------- Routines to set solver parameters ---------- */ 1403 1404 #undef __FUNCT__ 1405 #define __FUNCT__ "TSGetDuration" 1406 /*@ 1407 TSGetDuration - Gets the maximum number of timesteps to use and 1408 maximum time for iteration. 1409 1410 Not Collective 1411 1412 Input Parameters: 1413 + ts - the TS context obtained from TSCreate() 1414 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1415 - maxtime - final time to iterate to, or PETSC_NULL 1416 1417 Level: intermediate 1418 1419 .keywords: TS, timestep, get, maximum, iterations, time 1420 @*/ 1421 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1422 { 1423 PetscFunctionBegin; 1424 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1425 if (maxsteps) { 1426 PetscValidIntPointer(maxsteps,2); 1427 *maxsteps = ts->max_steps; 1428 } 1429 if (maxtime ) { 1430 PetscValidScalarPointer(maxtime,3); 1431 *maxtime = ts->max_time; 1432 } 1433 PetscFunctionReturn(0); 1434 } 1435 1436 #undef __FUNCT__ 1437 #define __FUNCT__ "TSSetDuration" 1438 /*@ 1439 TSSetDuration - Sets the maximum number of timesteps to use and 1440 maximum time for iteration. 1441 1442 Logically Collective on TS 1443 1444 Input Parameters: 1445 + ts - the TS context obtained from TSCreate() 1446 . maxsteps - maximum number of iterations to use 1447 - maxtime - final time to iterate to 1448 1449 Options Database Keys: 1450 . -ts_max_steps <maxsteps> - Sets maxsteps 1451 . -ts_max_time <maxtime> - Sets maxtime 1452 1453 Notes: 1454 The default maximum number of iterations is 5000. Default time is 5.0 1455 1456 Level: intermediate 1457 1458 .keywords: TS, timestep, set, maximum, iterations 1459 @*/ 1460 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1461 { 1462 PetscFunctionBegin; 1463 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1464 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1465 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1466 ts->max_steps = maxsteps; 1467 ts->max_time = maxtime; 1468 PetscFunctionReturn(0); 1469 } 1470 1471 #undef __FUNCT__ 1472 #define __FUNCT__ "TSSetSolution" 1473 /*@ 1474 TSSetSolution - Sets the initial solution vector 1475 for use by the TS routines. 1476 1477 Logically Collective on TS and Vec 1478 1479 Input Parameters: 1480 + ts - the TS context obtained from TSCreate() 1481 - x - the solution vector 1482 1483 Level: beginner 1484 1485 .keywords: TS, timestep, set, solution, initial conditions 1486 @*/ 1487 PetscErrorCode TSSetSolution(TS ts,Vec x) 1488 { 1489 PetscErrorCode ierr; 1490 1491 PetscFunctionBegin; 1492 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1493 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1494 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1495 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1496 ts->vec_sol = x; 1497 PetscFunctionReturn(0); 1498 } 1499 1500 #undef __FUNCT__ 1501 #define __FUNCT__ "TSSetPreStep" 1502 /*@C 1503 TSSetPreStep - Sets the general-purpose function 1504 called once at the beginning of each time step. 1505 1506 Logically Collective on TS 1507 1508 Input Parameters: 1509 + ts - The TS context obtained from TSCreate() 1510 - func - The function 1511 1512 Calling sequence of func: 1513 . func (TS ts); 1514 1515 Level: intermediate 1516 1517 .keywords: TS, timestep 1518 @*/ 1519 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1520 { 1521 PetscFunctionBegin; 1522 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1523 ts->ops->prestep = func; 1524 PetscFunctionReturn(0); 1525 } 1526 1527 #undef __FUNCT__ 1528 #define __FUNCT__ "TSPreStep" 1529 /*@C 1530 TSPreStep - Runs the user-defined pre-step function. 1531 1532 Collective on TS 1533 1534 Input Parameters: 1535 . ts - The TS context obtained from TSCreate() 1536 1537 Notes: 1538 TSPreStep() is typically used within time stepping implementations, 1539 so most users would not generally call this routine themselves. 1540 1541 Level: developer 1542 1543 .keywords: TS, timestep 1544 @*/ 1545 PetscErrorCode TSPreStep(TS ts) 1546 { 1547 PetscErrorCode ierr; 1548 1549 PetscFunctionBegin; 1550 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1551 if (ts->ops->prestep) { 1552 PetscStackPush("TS PreStep function"); 1553 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1554 PetscStackPop; 1555 } 1556 PetscFunctionReturn(0); 1557 } 1558 1559 #undef __FUNCT__ 1560 #define __FUNCT__ "TSSetPostStep" 1561 /*@C 1562 TSSetPostStep - Sets the general-purpose function 1563 called once at the end of each time step. 1564 1565 Logically Collective on TS 1566 1567 Input Parameters: 1568 + ts - The TS context obtained from TSCreate() 1569 - func - The function 1570 1571 Calling sequence of func: 1572 . func (TS ts); 1573 1574 Level: intermediate 1575 1576 .keywords: TS, timestep 1577 @*/ 1578 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1579 { 1580 PetscFunctionBegin; 1581 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1582 ts->ops->poststep = func; 1583 PetscFunctionReturn(0); 1584 } 1585 1586 #undef __FUNCT__ 1587 #define __FUNCT__ "TSPostStep" 1588 /*@C 1589 TSPostStep - Runs the user-defined post-step function. 1590 1591 Collective on TS 1592 1593 Input Parameters: 1594 . ts - The TS context obtained from TSCreate() 1595 1596 Notes: 1597 TSPostStep() is typically used within time stepping implementations, 1598 so most users would not generally call this routine themselves. 1599 1600 Level: developer 1601 1602 .keywords: TS, timestep 1603 @*/ 1604 PetscErrorCode TSPostStep(TS ts) 1605 { 1606 PetscErrorCode ierr; 1607 1608 PetscFunctionBegin; 1609 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1610 if (ts->ops->poststep) { 1611 PetscStackPush("TS PostStep function"); 1612 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1613 PetscStackPop; 1614 } 1615 PetscFunctionReturn(0); 1616 } 1617 1618 /* ------------ Routines to set performance monitoring options ----------- */ 1619 1620 #undef __FUNCT__ 1621 #define __FUNCT__ "TSMonitorSet" 1622 /*@C 1623 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1624 timestep to display the iteration's progress. 1625 1626 Logically Collective on TS 1627 1628 Input Parameters: 1629 + ts - the TS context obtained from TSCreate() 1630 . func - monitoring routine 1631 . mctx - [optional] user-defined context for private data for the 1632 monitor routine (use PETSC_NULL if no context is desired) 1633 - monitordestroy - [optional] routine that frees monitor context 1634 (may be PETSC_NULL) 1635 1636 Calling sequence of func: 1637 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1638 1639 + ts - the TS context 1640 . steps - iteration number 1641 . time - current time 1642 . x - current iterate 1643 - mctx - [optional] monitoring context 1644 1645 Notes: 1646 This routine adds an additional monitor to the list of monitors that 1647 already has been loaded. 1648 1649 Fortran notes: Only a single monitor function can be set for each TS object 1650 1651 Level: intermediate 1652 1653 .keywords: TS, timestep, set, monitor 1654 1655 .seealso: TSMonitorDefault(), TSMonitorCancel() 1656 @*/ 1657 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1658 { 1659 PetscFunctionBegin; 1660 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1661 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1662 ts->monitor[ts->numbermonitors] = monitor; 1663 ts->mdestroy[ts->numbermonitors] = mdestroy; 1664 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1665 PetscFunctionReturn(0); 1666 } 1667 1668 #undef __FUNCT__ 1669 #define __FUNCT__ "TSMonitorCancel" 1670 /*@C 1671 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1672 1673 Logically Collective on TS 1674 1675 Input Parameters: 1676 . ts - the TS context obtained from TSCreate() 1677 1678 Notes: 1679 There is no way to remove a single, specific monitor. 1680 1681 Level: intermediate 1682 1683 .keywords: TS, timestep, set, monitor 1684 1685 .seealso: TSMonitorDefault(), TSMonitorSet() 1686 @*/ 1687 PetscErrorCode TSMonitorCancel(TS ts) 1688 { 1689 PetscErrorCode ierr; 1690 PetscInt i; 1691 1692 PetscFunctionBegin; 1693 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1694 for (i=0; i<ts->numbermonitors; i++) { 1695 if (ts->mdestroy[i]) { 1696 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1697 } 1698 } 1699 ts->numbermonitors = 0; 1700 PetscFunctionReturn(0); 1701 } 1702 1703 #undef __FUNCT__ 1704 #define __FUNCT__ "TSMonitorDefault" 1705 /*@ 1706 TSMonitorDefault - Sets the Default monitor 1707 1708 Level: intermediate 1709 1710 .keywords: TS, set, monitor 1711 1712 .seealso: TSMonitorDefault(), TSMonitorSet() 1713 @*/ 1714 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1715 { 1716 PetscErrorCode ierr; 1717 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1718 1719 PetscFunctionBegin; 1720 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1721 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1722 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1723 PetscFunctionReturn(0); 1724 } 1725 1726 #undef __FUNCT__ 1727 #define __FUNCT__ "TSStep" 1728 /*@ 1729 TSStep - Steps the requested number of timesteps. 1730 1731 Collective on TS 1732 1733 Input Parameter: 1734 . ts - the TS context obtained from TSCreate() 1735 1736 Output Parameters: 1737 + steps - number of iterations until termination 1738 - ptime - time until termination 1739 1740 Level: beginner 1741 1742 .keywords: TS, timestep, solve 1743 1744 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1745 @*/ 1746 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1747 { 1748 PetscErrorCode ierr; 1749 1750 PetscFunctionBegin; 1751 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1752 1753 ierr = TSSetUp(ts);CHKERRQ(ierr); 1754 1755 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1756 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1757 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1758 1759 if (!PetscPreLoadingOn) { 1760 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1761 } 1762 PetscFunctionReturn(0); 1763 } 1764 1765 #undef __FUNCT__ 1766 #define __FUNCT__ "TSSolve" 1767 /*@ 1768 TSSolve - Steps the requested number of timesteps. 1769 1770 Collective on TS 1771 1772 Input Parameter: 1773 + ts - the TS context obtained from TSCreate() 1774 - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 1775 1776 Level: beginner 1777 1778 .keywords: TS, timestep, solve 1779 1780 .seealso: TSCreate(), TSSetSolution(), TSStep() 1781 @*/ 1782 PetscErrorCode TSSolve(TS ts, Vec x) 1783 { 1784 PetscInt steps; 1785 PetscReal ptime; 1786 PetscErrorCode ierr; 1787 1788 PetscFunctionBegin; 1789 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1790 /* set solution vector if provided */ 1791 if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 1792 /* reset time step and iteration counters */ 1793 ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 1794 /* steps the requested number of timesteps. */ 1795 ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 1796 PetscFunctionReturn(0); 1797 } 1798 1799 #undef __FUNCT__ 1800 #define __FUNCT__ "TSMonitor" 1801 /* 1802 Runs the user provided monitor routines, if they exists. 1803 */ 1804 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1805 { 1806 PetscErrorCode ierr; 1807 PetscInt i,n = ts->numbermonitors; 1808 1809 PetscFunctionBegin; 1810 for (i=0; i<n; i++) { 1811 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1812 } 1813 PetscFunctionReturn(0); 1814 } 1815 1816 /* ------------------------------------------------------------------------*/ 1817 1818 #undef __FUNCT__ 1819 #define __FUNCT__ "TSMonitorLGCreate" 1820 /*@C 1821 TSMonitorLGCreate - Creates a line graph context for use with 1822 TS to monitor convergence of preconditioned residual norms. 1823 1824 Collective on TS 1825 1826 Input Parameters: 1827 + host - the X display to open, or null for the local machine 1828 . label - the title to put in the title bar 1829 . x, y - the screen coordinates of the upper left coordinate of the window 1830 - m, n - the screen width and height in pixels 1831 1832 Output Parameter: 1833 . draw - the drawing context 1834 1835 Options Database Key: 1836 . -ts_monitor_draw - automatically sets line graph monitor 1837 1838 Notes: 1839 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1840 1841 Level: intermediate 1842 1843 .keywords: TS, monitor, line graph, residual, seealso 1844 1845 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1846 1847 @*/ 1848 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1849 { 1850 PetscDraw win; 1851 PetscErrorCode ierr; 1852 1853 PetscFunctionBegin; 1854 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1855 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1856 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1857 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1858 1859 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1860 PetscFunctionReturn(0); 1861 } 1862 1863 #undef __FUNCT__ 1864 #define __FUNCT__ "TSMonitorLG" 1865 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1866 { 1867 PetscDrawLG lg = (PetscDrawLG) monctx; 1868 PetscReal x,y = ptime; 1869 PetscErrorCode ierr; 1870 1871 PetscFunctionBegin; 1872 if (!monctx) { 1873 MPI_Comm comm; 1874 PetscViewer viewer; 1875 1876 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1877 viewer = PETSC_VIEWER_DRAW_(comm); 1878 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1879 } 1880 1881 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1882 x = (PetscReal)n; 1883 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1884 if (n < 20 || (n % 5)) { 1885 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1886 } 1887 PetscFunctionReturn(0); 1888 } 1889 1890 #undef __FUNCT__ 1891 #define __FUNCT__ "TSMonitorLGDestroy" 1892 /*@C 1893 TSMonitorLGDestroy - Destroys a line graph context that was created 1894 with TSMonitorLGCreate(). 1895 1896 Collective on PetscDrawLG 1897 1898 Input Parameter: 1899 . draw - the drawing context 1900 1901 Level: intermediate 1902 1903 .keywords: TS, monitor, line graph, destroy 1904 1905 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1906 @*/ 1907 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1908 { 1909 PetscDraw draw; 1910 PetscErrorCode ierr; 1911 1912 PetscFunctionBegin; 1913 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 1914 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1915 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1916 PetscFunctionReturn(0); 1917 } 1918 1919 #undef __FUNCT__ 1920 #define __FUNCT__ "TSGetTime" 1921 /*@ 1922 TSGetTime - Gets the current time. 1923 1924 Not Collective 1925 1926 Input Parameter: 1927 . ts - the TS context obtained from TSCreate() 1928 1929 Output Parameter: 1930 . t - the current time 1931 1932 Level: beginner 1933 1934 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1935 1936 .keywords: TS, get, time 1937 @*/ 1938 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1939 { 1940 PetscFunctionBegin; 1941 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1942 PetscValidDoublePointer(t,2); 1943 *t = ts->ptime; 1944 PetscFunctionReturn(0); 1945 } 1946 1947 #undef __FUNCT__ 1948 #define __FUNCT__ "TSSetTime" 1949 /*@ 1950 TSSetTime - Allows one to reset the time. 1951 1952 Logically Collective on TS 1953 1954 Input Parameters: 1955 + ts - the TS context obtained from TSCreate() 1956 - time - the time 1957 1958 Level: intermediate 1959 1960 .seealso: TSGetTime(), TSSetDuration() 1961 1962 .keywords: TS, set, time 1963 @*/ 1964 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1965 { 1966 PetscFunctionBegin; 1967 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1968 PetscValidLogicalCollectiveReal(ts,t,2); 1969 ts->ptime = t; 1970 PetscFunctionReturn(0); 1971 } 1972 1973 #undef __FUNCT__ 1974 #define __FUNCT__ "TSSetOptionsPrefix" 1975 /*@C 1976 TSSetOptionsPrefix - Sets the prefix used for searching for all 1977 TS options in the database. 1978 1979 Logically Collective on TS 1980 1981 Input Parameter: 1982 + ts - The TS context 1983 - prefix - The prefix to prepend to all option names 1984 1985 Notes: 1986 A hyphen (-) must NOT be given at the beginning of the prefix name. 1987 The first character of all runtime options is AUTOMATICALLY the 1988 hyphen. 1989 1990 Level: advanced 1991 1992 .keywords: TS, set, options, prefix, database 1993 1994 .seealso: TSSetFromOptions() 1995 1996 @*/ 1997 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1998 { 1999 PetscErrorCode ierr; 2000 SNES snes; 2001 2002 PetscFunctionBegin; 2003 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2004 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2005 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2006 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2007 PetscFunctionReturn(0); 2008 } 2009 2010 2011 #undef __FUNCT__ 2012 #define __FUNCT__ "TSAppendOptionsPrefix" 2013 /*@C 2014 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 2015 TS options in the database. 2016 2017 Logically Collective on TS 2018 2019 Input Parameter: 2020 + ts - The TS context 2021 - prefix - The prefix to prepend to all option names 2022 2023 Notes: 2024 A hyphen (-) must NOT be given at the beginning of the prefix name. 2025 The first character of all runtime options is AUTOMATICALLY the 2026 hyphen. 2027 2028 Level: advanced 2029 2030 .keywords: TS, append, options, prefix, database 2031 2032 .seealso: TSGetOptionsPrefix() 2033 2034 @*/ 2035 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2036 { 2037 PetscErrorCode ierr; 2038 SNES snes; 2039 2040 PetscFunctionBegin; 2041 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2042 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2043 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2044 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2045 PetscFunctionReturn(0); 2046 } 2047 2048 #undef __FUNCT__ 2049 #define __FUNCT__ "TSGetOptionsPrefix" 2050 /*@C 2051 TSGetOptionsPrefix - Sets the prefix used for searching for all 2052 TS options in the database. 2053 2054 Not Collective 2055 2056 Input Parameter: 2057 . ts - The TS context 2058 2059 Output Parameter: 2060 . prefix - A pointer to the prefix string used 2061 2062 Notes: On the fortran side, the user should pass in a string 'prifix' of 2063 sufficient length to hold the prefix. 2064 2065 Level: intermediate 2066 2067 .keywords: TS, get, options, prefix, database 2068 2069 .seealso: TSAppendOptionsPrefix() 2070 @*/ 2071 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2072 { 2073 PetscErrorCode ierr; 2074 2075 PetscFunctionBegin; 2076 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2077 PetscValidPointer(prefix,2); 2078 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2079 PetscFunctionReturn(0); 2080 } 2081 2082 #undef __FUNCT__ 2083 #define __FUNCT__ "TSGetRHSJacobian" 2084 /*@C 2085 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2086 2087 Not Collective, but parallel objects are returned if TS is parallel 2088 2089 Input Parameter: 2090 . ts - The TS context obtained from TSCreate() 2091 2092 Output Parameters: 2093 + J - The Jacobian J of F, where U_t = F(U,t) 2094 . M - The preconditioner matrix, usually the same as J 2095 . func - Function to compute the Jacobian of the RHS 2096 - ctx - User-defined context for Jacobian evaluation routine 2097 2098 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2099 2100 Level: intermediate 2101 2102 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2103 2104 .keywords: TS, timestep, get, matrix, Jacobian 2105 @*/ 2106 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2107 { 2108 PetscErrorCode ierr; 2109 SNES snes; 2110 2111 PetscFunctionBegin; 2112 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2113 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2114 if (func) *func = ts->ops->rhsjacobian; 2115 if (ctx) *ctx = ts->jacP; 2116 PetscFunctionReturn(0); 2117 } 2118 2119 #undef __FUNCT__ 2120 #define __FUNCT__ "TSGetIJacobian" 2121 /*@C 2122 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2123 2124 Not Collective, but parallel objects are returned if TS is parallel 2125 2126 Input Parameter: 2127 . ts - The TS context obtained from TSCreate() 2128 2129 Output Parameters: 2130 + A - The Jacobian of F(t,U,U_t) 2131 . B - The preconditioner matrix, often the same as A 2132 . f - The function to compute the matrices 2133 - ctx - User-defined context for Jacobian evaluation routine 2134 2135 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2136 2137 Level: advanced 2138 2139 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2140 2141 .keywords: TS, timestep, get, matrix, Jacobian 2142 @*/ 2143 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2144 { 2145 PetscErrorCode ierr; 2146 SNES snes; 2147 2148 PetscFunctionBegin; 2149 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2150 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2151 if (f) *f = ts->ops->ijacobian; 2152 if (ctx) *ctx = ts->jacP; 2153 PetscFunctionReturn(0); 2154 } 2155 2156 #undef __FUNCT__ 2157 #define __FUNCT__ "TSMonitorSolution" 2158 /*@C 2159 TSMonitorSolution - Monitors progress of the TS solvers by calling 2160 VecView() for the solution at each timestep 2161 2162 Collective on TS 2163 2164 Input Parameters: 2165 + ts - the TS context 2166 . step - current time-step 2167 . ptime - current time 2168 - dummy - either a viewer or PETSC_NULL 2169 2170 Level: intermediate 2171 2172 .keywords: TS, vector, monitor, view 2173 2174 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2175 @*/ 2176 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2177 { 2178 PetscErrorCode ierr; 2179 PetscViewer viewer = (PetscViewer) dummy; 2180 2181 PetscFunctionBegin; 2182 if (!dummy) { 2183 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2184 } 2185 ierr = VecView(x,viewer);CHKERRQ(ierr); 2186 PetscFunctionReturn(0); 2187 } 2188 2189 2190 #undef __FUNCT__ 2191 #define __FUNCT__ "TSSetDM" 2192 /*@ 2193 TSSetDM - Sets the DM that may be used by some preconditioners 2194 2195 Logically Collective on TS and DM 2196 2197 Input Parameters: 2198 + ts - the preconditioner context 2199 - dm - the dm 2200 2201 Level: intermediate 2202 2203 2204 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2205 @*/ 2206 PetscErrorCode TSSetDM(TS ts,DM dm) 2207 { 2208 PetscErrorCode ierr; 2209 SNES snes; 2210 2211 PetscFunctionBegin; 2212 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2213 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2214 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2215 ts->dm = dm; 2216 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2217 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2218 PetscFunctionReturn(0); 2219 } 2220 2221 #undef __FUNCT__ 2222 #define __FUNCT__ "TSGetDM" 2223 /*@ 2224 TSGetDM - Gets the DM that may be used by some preconditioners 2225 2226 Not Collective 2227 2228 Input Parameter: 2229 . ts - the preconditioner context 2230 2231 Output Parameter: 2232 . dm - the dm 2233 2234 Level: intermediate 2235 2236 2237 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2238 @*/ 2239 PetscErrorCode TSGetDM(TS ts,DM *dm) 2240 { 2241 PetscFunctionBegin; 2242 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2243 *dm = ts->dm; 2244 PetscFunctionReturn(0); 2245 } 2246 2247 #undef __FUNCT__ 2248 #define __FUNCT__ "SNESTSFormFunction" 2249 /*@ 2250 SNESTSFormFunction - Function to evaluate nonlinear residual 2251 2252 Logically Collective on SNES 2253 2254 Input Parameter: 2255 + snes - nonlinear solver 2256 . X - the current state at which to evaluate the residual 2257 - ctx - user context, must be a TS 2258 2259 Output Parameter: 2260 . F - the nonlinear residual 2261 2262 Notes: 2263 This function is not normally called by users and is automatically registered with the SNES used by TS. 2264 It is most frequently passed to MatFDColoringSetFunction(). 2265 2266 Level: advanced 2267 2268 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2269 @*/ 2270 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2271 { 2272 TS ts = (TS)ctx; 2273 PetscErrorCode ierr; 2274 2275 PetscFunctionBegin; 2276 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2277 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2278 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2279 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2280 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2281 PetscFunctionReturn(0); 2282 } 2283 2284 #undef __FUNCT__ 2285 #define __FUNCT__ "SNESTSFormJacobian" 2286 /*@ 2287 SNESTSFormJacobian - Function to evaluate the Jacobian 2288 2289 Collective on SNES 2290 2291 Input Parameter: 2292 + snes - nonlinear solver 2293 . X - the current state at which to evaluate the residual 2294 - ctx - user context, must be a TS 2295 2296 Output Parameter: 2297 + A - the Jacobian 2298 . B - the preconditioning matrix (may be the same as A) 2299 - flag - indicates any structure change in the matrix 2300 2301 Notes: 2302 This function is not normally called by users and is automatically registered with the SNES used by TS. 2303 2304 Level: developer 2305 2306 .seealso: SNESSetJacobian() 2307 @*/ 2308 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2309 { 2310 TS ts = (TS)ctx; 2311 PetscErrorCode ierr; 2312 2313 PetscFunctionBegin; 2314 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2315 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2316 PetscValidPointer(A,3); 2317 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2318 PetscValidPointer(B,4); 2319 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2320 PetscValidPointer(flag,5); 2321 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2322 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2323 PetscFunctionReturn(0); 2324 } 2325 2326 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2327 #include <mex.h> 2328 2329 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2330 2331 #undef __FUNCT__ 2332 #define __FUNCT__ "TSComputeFunction_Matlab" 2333 /* 2334 TSComputeFunction_Matlab - Calls the function that has been set with 2335 TSSetFunctionMatlab(). 2336 2337 Collective on TS 2338 2339 Input Parameters: 2340 + snes - the TS context 2341 - x - input vector 2342 2343 Output Parameter: 2344 . y - function vector, as set by TSSetFunction() 2345 2346 Notes: 2347 TSComputeFunction() is typically used within nonlinear solvers 2348 implementations, so most users would not generally call this routine 2349 themselves. 2350 2351 Level: developer 2352 2353 .keywords: TS, nonlinear, compute, function 2354 2355 .seealso: TSSetFunction(), TSGetFunction() 2356 */ 2357 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2358 { 2359 PetscErrorCode ierr; 2360 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2361 int nlhs = 1,nrhs = 7; 2362 mxArray *plhs[1],*prhs[7]; 2363 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2364 2365 PetscFunctionBegin; 2366 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2367 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2368 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2369 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2370 PetscCheckSameComm(snes,1,x,3); 2371 PetscCheckSameComm(snes,1,y,5); 2372 2373 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2374 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2375 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2376 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2377 prhs[0] = mxCreateDoubleScalar((double)ls); 2378 prhs[1] = mxCreateDoubleScalar(time); 2379 prhs[2] = mxCreateDoubleScalar((double)lx); 2380 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2381 prhs[4] = mxCreateDoubleScalar((double)ly); 2382 prhs[5] = mxCreateString(sctx->funcname); 2383 prhs[6] = sctx->ctx; 2384 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2385 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2386 mxDestroyArray(prhs[0]); 2387 mxDestroyArray(prhs[1]); 2388 mxDestroyArray(prhs[2]); 2389 mxDestroyArray(prhs[3]); 2390 mxDestroyArray(prhs[4]); 2391 mxDestroyArray(prhs[5]); 2392 mxDestroyArray(plhs[0]); 2393 PetscFunctionReturn(0); 2394 } 2395 2396 2397 #undef __FUNCT__ 2398 #define __FUNCT__ "TSSetFunctionMatlab" 2399 /* 2400 TSSetFunctionMatlab - Sets the function evaluation routine and function 2401 vector for use by the TS routines in solving ODEs 2402 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2403 2404 Logically Collective on TS 2405 2406 Input Parameters: 2407 + ts - the TS context 2408 - func - function evaluation routine 2409 2410 Calling sequence of func: 2411 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2412 2413 Level: beginner 2414 2415 .keywords: TS, nonlinear, set, function 2416 2417 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2418 */ 2419 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2420 { 2421 PetscErrorCode ierr; 2422 TSMatlabContext *sctx; 2423 2424 PetscFunctionBegin; 2425 /* currently sctx is memory bleed */ 2426 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2427 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2428 /* 2429 This should work, but it doesn't 2430 sctx->ctx = ctx; 2431 mexMakeArrayPersistent(sctx->ctx); 2432 */ 2433 sctx->ctx = mxDuplicateArray(ctx); 2434 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2435 PetscFunctionReturn(0); 2436 } 2437 2438 #undef __FUNCT__ 2439 #define __FUNCT__ "TSComputeJacobian_Matlab" 2440 /* 2441 TSComputeJacobian_Matlab - Calls the function that has been set with 2442 TSSetJacobianMatlab(). 2443 2444 Collective on TS 2445 2446 Input Parameters: 2447 + snes - the TS context 2448 . x - input vector 2449 . A, B - the matrices 2450 - ctx - user context 2451 2452 Output Parameter: 2453 . flag - structure of the matrix 2454 2455 Level: developer 2456 2457 .keywords: TS, nonlinear, compute, function 2458 2459 .seealso: TSSetFunction(), TSGetFunction() 2460 @*/ 2461 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2462 { 2463 PetscErrorCode ierr; 2464 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2465 int nlhs = 2,nrhs = 9; 2466 mxArray *plhs[2],*prhs[9]; 2467 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2468 2469 PetscFunctionBegin; 2470 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2471 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2472 2473 /* call Matlab function in ctx with arguments x and y */ 2474 2475 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2476 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2477 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2478 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2479 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2480 prhs[0] = mxCreateDoubleScalar((double)ls); 2481 prhs[1] = mxCreateDoubleScalar((double)time); 2482 prhs[2] = mxCreateDoubleScalar((double)lx); 2483 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2484 prhs[4] = mxCreateDoubleScalar((double)shift); 2485 prhs[5] = mxCreateDoubleScalar((double)lA); 2486 prhs[6] = mxCreateDoubleScalar((double)lB); 2487 prhs[7] = mxCreateString(sctx->funcname); 2488 prhs[8] = sctx->ctx; 2489 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2490 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2491 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2492 mxDestroyArray(prhs[0]); 2493 mxDestroyArray(prhs[1]); 2494 mxDestroyArray(prhs[2]); 2495 mxDestroyArray(prhs[3]); 2496 mxDestroyArray(prhs[4]); 2497 mxDestroyArray(prhs[5]); 2498 mxDestroyArray(prhs[6]); 2499 mxDestroyArray(prhs[7]); 2500 mxDestroyArray(plhs[0]); 2501 mxDestroyArray(plhs[1]); 2502 PetscFunctionReturn(0); 2503 } 2504 2505 2506 #undef __FUNCT__ 2507 #define __FUNCT__ "TSSetJacobianMatlab" 2508 /* 2509 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2510 vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function 2511 2512 Logically Collective on TS 2513 2514 Input Parameters: 2515 + snes - the TS context 2516 . A,B - Jacobian matrices 2517 . func - function evaluation routine 2518 - ctx - user context 2519 2520 Calling sequence of func: 2521 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2522 2523 2524 Level: developer 2525 2526 .keywords: TS, nonlinear, set, function 2527 2528 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2529 */ 2530 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2531 { 2532 PetscErrorCode ierr; 2533 TSMatlabContext *sctx; 2534 2535 PetscFunctionBegin; 2536 /* currently sctx is memory bleed */ 2537 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2538 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2539 /* 2540 This should work, but it doesn't 2541 sctx->ctx = ctx; 2542 mexMakeArrayPersistent(sctx->ctx); 2543 */ 2544 sctx->ctx = mxDuplicateArray(ctx); 2545 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2546 PetscFunctionReturn(0); 2547 } 2548 2549 #undef __FUNCT__ 2550 #define __FUNCT__ "TSMonitor_Matlab" 2551 /* 2552 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2553 2554 Collective on TS 2555 2556 .seealso: TSSetFunction(), TSGetFunction() 2557 @*/ 2558 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2559 { 2560 PetscErrorCode ierr; 2561 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2562 int nlhs = 1,nrhs = 6; 2563 mxArray *plhs[1],*prhs[6]; 2564 long long int lx = 0,ls = 0; 2565 2566 PetscFunctionBegin; 2567 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2568 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2569 2570 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2571 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2572 prhs[0] = mxCreateDoubleScalar((double)ls); 2573 prhs[1] = mxCreateDoubleScalar((double)it); 2574 prhs[2] = mxCreateDoubleScalar((double)time); 2575 prhs[3] = mxCreateDoubleScalar((double)lx); 2576 prhs[4] = mxCreateString(sctx->funcname); 2577 prhs[5] = sctx->ctx; 2578 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2579 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2580 mxDestroyArray(prhs[0]); 2581 mxDestroyArray(prhs[1]); 2582 mxDestroyArray(prhs[2]); 2583 mxDestroyArray(prhs[3]); 2584 mxDestroyArray(prhs[4]); 2585 mxDestroyArray(plhs[0]); 2586 PetscFunctionReturn(0); 2587 } 2588 2589 2590 #undef __FUNCT__ 2591 #define __FUNCT__ "TSMonitorSetMatlab" 2592 /* 2593 TSMonitorSetMatlab - Sets the monitor function from Matlab 2594 2595 Level: developer 2596 2597 .keywords: TS, nonlinear, set, function 2598 2599 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2600 */ 2601 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2602 { 2603 PetscErrorCode ierr; 2604 TSMatlabContext *sctx; 2605 2606 PetscFunctionBegin; 2607 /* currently sctx is memory bleed */ 2608 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2609 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2610 /* 2611 This should work, but it doesn't 2612 sctx->ctx = ctx; 2613 mexMakeArrayPersistent(sctx->ctx); 2614 */ 2615 sctx->ctx = mxDuplicateArray(ctx); 2616 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2617 PetscFunctionReturn(0); 2618 } 2619 2620 #endif 2621