1 #include "src/ts/tsimpl.h" /*I "petscts.h" I*/ 2 3 /* Logging support */ 4 PetscCookie TS_COOKIE = 0; 5 PetscEvent TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0; 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "TSSetTypeFromOptions" 9 /* 10 TSSetTypeFromOptions - Sets the type of ts from user options. 11 12 Collective on TS 13 14 Input Parameter: 15 . ts - The ts 16 17 Level: intermediate 18 19 .keywords: TS, set, options, database, type 20 .seealso: TSSetFromOptions(), TSSetType() 21 */ 22 static PetscErrorCode TSSetTypeFromOptions(TS ts) 23 { 24 PetscTruth opt; 25 const char *defaultType; 26 char typeName[256]; 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 if (ts->type_name != PETSC_NULL) { 31 defaultType = ts->type_name; 32 } else { 33 defaultType = TS_EULER; 34 } 35 36 if (!TSRegisterAllCalled) { 37 ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr); 38 } 39 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 40 if (opt) { 41 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 42 } else { 43 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 44 } 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "TSSetFromOptions" 50 /*@ 51 TSSetFromOptions - Sets various TS parameters from user options. 52 53 Collective on TS 54 55 Input Parameter: 56 . ts - the TS context obtained from TSCreate() 57 58 Options Database Keys: 59 + -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON 60 . -ts_max_steps maxsteps - maximum number of time-steps to take 61 . -ts_max_time time - maximum time to compute to 62 . -ts_dt dt - initial time step 63 . -ts_monitor - print information at each timestep 64 - -ts_xmonitor - plot information at each timestep 65 66 Level: beginner 67 68 .keywords: TS, timestep, set, options, database 69 70 .seealso: TSGetType 71 @*/ 72 PetscErrorCode TSSetFromOptions(TS ts) 73 { 74 PetscReal dt; 75 PetscTruth opt; 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 80 ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");CHKERRQ(ierr); 81 82 /* Handle generic TS options */ 83 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 84 ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 85 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 87 if (opt) { 88 ts->initial_time_step = ts->time_step = dt; 89 } 90 91 /* Monitor options */ 92 ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);CHKERRQ(ierr); 93 if (opt) { 94 ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 95 } 96 ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);CHKERRQ(ierr); 97 if (opt) { 98 ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 99 } 100 ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);CHKERRQ(ierr); 101 if (opt) { 102 ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 103 } 104 105 /* Handle TS type options */ 106 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 107 108 /* Handle specific TS options */ 109 if (ts->ops->setfromoptions != PETSC_NULL) { 110 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 111 } 112 ierr = PetscOptionsEnd();CHKERRQ(ierr); 113 114 /* Handle subobject options */ 115 switch(ts->problem_type) { 116 /* Should check for implicit/explicit */ 117 case TS_LINEAR: 118 if (ts->ksp != PETSC_NULL) { 119 ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr); 120 } 121 break; 122 case TS_NONLINEAR: 123 if (ts->snes != PETSC_NULL) { 124 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 125 } 126 break; 127 default: 128 SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type); 129 } 130 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "TSViewFromOptions" 136 /*@ 137 TSViewFromOptions - This function visualizes the ts based upon user options. 138 139 Collective on TS 140 141 Input Parameter: 142 . ts - The ts 143 144 Level: intermediate 145 146 .keywords: TS, view, options, database 147 .seealso: TSSetFromOptions(), TSView() 148 @*/ 149 PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 150 { 151 PetscViewer viewer; 152 PetscDraw draw; 153 PetscTruth opt; 154 char typeName[1024]; 155 char fileName[PETSC_MAX_PATH_LEN]; 156 size_t len; 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);CHKERRQ(ierr); 161 if (opt) { 162 ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);CHKERRQ(ierr); 163 ierr = PetscStrlen(typeName, &len);CHKERRQ(ierr); 164 if (len > 0) { 165 ierr = PetscViewerCreate(ts->comm, &viewer);CHKERRQ(ierr); 166 ierr = PetscViewerSetType(viewer, typeName);CHKERRQ(ierr); 167 ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);CHKERRQ(ierr); 168 if (opt) { 169 ierr = PetscViewerSetFilename(viewer, fileName);CHKERRQ(ierr); 170 } else { 171 ierr = PetscViewerSetFilename(viewer, ts->name);CHKERRQ(ierr); 172 } 173 ierr = TSView(ts, viewer);CHKERRQ(ierr); 174 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 175 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 176 } else { 177 ierr = TSView(ts, PETSC_NULL);CHKERRQ(ierr); 178 } 179 } 180 ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr); 181 if (opt) { 182 ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 183 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 184 if (title) { 185 ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 186 } else { 187 ierr = PetscObjectName((PetscObject) ts);CHKERRQ(ierr); 188 ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr); 189 } 190 ierr = TSView(ts, viewer);CHKERRQ(ierr); 191 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 192 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 193 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 194 } 195 PetscFunctionReturn(0); 196 } 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "TSComputeRHSJacobian" 200 /*@ 201 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 202 set with TSSetRHSJacobian(). 203 204 Collective on TS and Vec 205 206 Input Parameters: 207 + ts - the SNES context 208 . t - current timestep 209 - x - input vector 210 211 Output Parameters: 212 + A - Jacobian matrix 213 . B - optional preconditioning matrix 214 - flag - flag indicating matrix structure 215 216 Notes: 217 Most users should not need to explicitly call this routine, as it 218 is used internally within the nonlinear solvers. 219 220 See KSPSetOperators() for important information about setting the 221 flag parameter. 222 223 TSComputeJacobian() is valid only for TS_NONLINEAR 224 225 Level: developer 226 227 .keywords: SNES, compute, Jacobian, matrix 228 229 .seealso: TSSetRHSJacobian(), KSPSetOperators() 230 @*/ 231 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 232 { 233 PetscErrorCode ierr; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 237 PetscValidHeaderSpecific(X,VEC_COOKIE,3); 238 PetscCheckSameComm(ts,1,X,3); 239 if (ts->problem_type != TS_NONLINEAR) { 240 SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only"); 241 } 242 if (ts->ops->rhsjacobian) { 243 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 244 *flg = DIFFERENT_NONZERO_PATTERN; 245 PetscStackPush("TS user Jacobian function"); 246 ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 247 PetscStackPop; 248 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 249 /* make sure user returned a correct Jacobian and preconditioner */ 250 PetscValidHeaderSpecific(*A,MAT_COOKIE,4); 251 PetscValidHeaderSpecific(*B,MAT_COOKIE,5); 252 } else { 253 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 254 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255 if (*A != *B) { 256 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 257 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258 } 259 } 260 PetscFunctionReturn(0); 261 } 262 263 #undef __FUNCT__ 264 #define __FUNCT__ "TSComputeRHSFunction" 265 /* 266 TSComputeRHSFunction - Evaluates the right-hand-side function. 267 268 Note: If the user did not provide a function but merely a matrix, 269 this routine applies the matrix. 270 */ 271 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 272 { 273 PetscErrorCode ierr; 274 275 PetscFunctionBegin; 276 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 277 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 278 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 279 280 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 281 if (ts->ops->rhsfunction) { 282 PetscStackPush("TS user right-hand-side function"); 283 ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 284 PetscStackPop; 285 } else { 286 if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 287 MatStructure flg; 288 PetscStackPush("TS user right-hand-side matrix function"); 289 ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 290 PetscStackPop; 291 } 292 ierr = MatMult(ts->A,x,y);CHKERRQ(ierr); 293 } 294 295 /* apply user-provided boundary conditions (only needed if these are time dependent) */ 296 ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr); 297 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 298 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "TSSetRHSFunction" 304 /*@C 305 TSSetRHSFunction - Sets the routine for evaluating the function, 306 F(t,u), where U_t = F(t,u). 307 308 Collective on TS 309 310 Input Parameters: 311 + ts - the TS context obtained from TSCreate() 312 . f - routine for evaluating the right-hand-side function 313 - ctx - [optional] user-defined context for private data for the 314 function evaluation routine (may be PETSC_NULL) 315 316 Calling sequence of func: 317 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 318 319 + t - current timestep 320 . u - input vector 321 . F - function vector 322 - ctx - [optional] user-defined function context 323 324 Important: 325 The user MUST call either this routine or TSSetRHSMatrix(). 326 327 Level: beginner 328 329 .keywords: TS, timestep, set, right-hand-side, function 330 331 .seealso: TSSetRHSMatrix() 332 @*/ 333 PetscErrorCode TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 334 { 335 PetscFunctionBegin; 336 337 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 338 if (ts->problem_type == TS_LINEAR) { 339 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 340 } 341 ts->ops->rhsfunction = f; 342 ts->funP = ctx; 343 PetscFunctionReturn(0); 344 } 345 346 #undef __FUNCT__ 347 #define __FUNCT__ "TSSetRHSMatrix" 348 /*@C 349 TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U. 350 Also sets the location to store A. 351 352 Collective on TS 353 354 Input Parameters: 355 + ts - the TS context obtained from TSCreate() 356 . A - matrix 357 . B - preconditioner matrix (usually same as A) 358 . f - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 359 if A is not a function of t. 360 - ctx - [optional] user-defined context for private data for the 361 matrix evaluation routine (may be PETSC_NULL) 362 363 Calling sequence of func: 364 $ func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 365 366 + t - current timestep 367 . A - matrix A, where U_t = A(t) U 368 . B - preconditioner matrix, usually the same as A 369 . flag - flag indicating information about the preconditioner matrix 370 structure (same as flag in KSPSetOperators()) 371 - ctx - [optional] user-defined context for matrix evaluation routine 372 373 Notes: 374 See KSPSetOperators() for important information about setting the flag 375 output parameter in the routine func(). Be sure to read this information! 376 377 The routine func() takes Mat * as the matrix arguments rather than Mat. 378 This allows the matrix evaluation routine to replace A and/or B with a 379 completely new new matrix structure (not just different matrix elements) 380 when appropriate, for instance, if the nonzero structure is changing 381 throughout the global iterations. 382 383 Important: 384 The user MUST call either this routine or TSSetRHSFunction(). 385 386 Level: beginner 387 388 .keywords: TS, timestep, set, right-hand-side, matrix 389 390 .seealso: TSSetRHSFunction() 391 @*/ 392 PetscErrorCode TSSetRHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx) 393 { 394 PetscFunctionBegin; 395 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 396 PetscValidHeaderSpecific(A,MAT_COOKIE,2); 397 PetscValidHeaderSpecific(B,MAT_COOKIE,3); 398 PetscCheckSameComm(ts,1,A,2); 399 PetscCheckSameComm(ts,1,B,3); 400 if (ts->problem_type == TS_NONLINEAR) { 401 SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()"); 402 } 403 404 ts->ops->rhsmatrix = f; 405 ts->jacP = ctx; 406 ts->A = A; 407 ts->B = B; 408 409 PetscFunctionReturn(0); 410 } 411 412 #undef __FUNCT__ 413 #define __FUNCT__ "TSSetRHSJacobian" 414 /*@C 415 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 416 where U_t = F(U,t), as well as the location to store the matrix. 417 Use TSSetRHSMatrix() for linear problems. 418 419 Collective on TS 420 421 Input Parameters: 422 + ts - the TS context obtained from TSCreate() 423 . A - Jacobian matrix 424 . B - preconditioner matrix (usually same as A) 425 . f - the Jacobian evaluation routine 426 - ctx - [optional] user-defined context for private data for the 427 Jacobian evaluation routine (may be PETSC_NULL) 428 429 Calling sequence of func: 430 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 431 432 + t - current timestep 433 . u - input vector 434 . A - matrix A, where U_t = A(t)u 435 . B - preconditioner matrix, usually the same as A 436 . flag - flag indicating information about the preconditioner matrix 437 structure (same as flag in KSPSetOperators()) 438 - ctx - [optional] user-defined context for matrix evaluation routine 439 440 Notes: 441 See KSPSetOperators() for important information about setting the flag 442 output parameter in the routine func(). Be sure to read this information! 443 444 The routine func() takes Mat * as the matrix arguments rather than Mat. 445 This allows the matrix evaluation routine to replace A and/or B with a 446 completely new new matrix structure (not just different matrix elements) 447 when appropriate, for instance, if the nonzero structure is changing 448 throughout the global iterations. 449 450 Level: beginner 451 452 .keywords: TS, timestep, set, right-hand-side, Jacobian 453 454 .seealso: TSDefaultComputeJacobianColor(), 455 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix() 456 457 @*/ 458 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 459 { 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 462 PetscValidHeaderSpecific(A,MAT_COOKIE,2); 463 PetscValidHeaderSpecific(B,MAT_COOKIE,3); 464 PetscCheckSameComm(ts,1,A,2); 465 PetscCheckSameComm(ts,1,B,3); 466 if (ts->problem_type != TS_NONLINEAR) { 467 SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()"); 468 } 469 470 ts->ops->rhsjacobian = f; 471 ts->jacP = ctx; 472 ts->A = A; 473 ts->B = B; 474 PetscFunctionReturn(0); 475 } 476 477 #undef __FUNCT__ 478 #define __FUNCT__ "TSComputeRHSBoundaryConditions" 479 /* 480 TSComputeRHSBoundaryConditions - Evaluates the boundary condition function. 481 482 Note: If the user did not provide a function but merely a matrix, 483 this routine applies the matrix. 484 */ 485 PetscErrorCode TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x) 486 { 487 PetscErrorCode ierr; 488 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 491 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 492 PetscCheckSameComm(ts,1,x,3); 493 494 if (ts->ops->rhsbc) { 495 PetscStackPush("TS user boundary condition function"); 496 ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr); 497 PetscStackPop; 498 PetscFunctionReturn(0); 499 } 500 501 PetscFunctionReturn(0); 502 } 503 504 #undef __FUNCT__ 505 #define __FUNCT__ "TSSetRHSBoundaryConditions" 506 /*@C 507 TSSetRHSBoundaryConditions - Sets the routine for evaluating the function, 508 boundary conditions for the function F. 509 510 Collective on TS 511 512 Input Parameters: 513 + ts - the TS context obtained from TSCreate() 514 . f - routine for evaluating the boundary condition function 515 - ctx - [optional] user-defined context for private data for the 516 function evaluation routine (may be PETSC_NULL) 517 518 Calling sequence of func: 519 $ func (TS ts,PetscReal t,Vec F,void *ctx); 520 521 + t - current timestep 522 . F - function vector 523 - ctx - [optional] user-defined function context 524 525 Level: intermediate 526 527 .keywords: TS, timestep, set, boundary conditions, function 528 @*/ 529 PetscErrorCode TSSetRHSBoundaryConditions(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx) 530 { 531 PetscFunctionBegin; 532 533 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 534 if (ts->problem_type != TS_LINEAR) { 535 SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only"); 536 } 537 ts->ops->rhsbc = f; 538 ts->bcP = ctx; 539 PetscFunctionReturn(0); 540 } 541 542 #undef __FUNCT__ 543 #define __FUNCT__ "TSView" 544 /*@C 545 TSView - Prints the TS data structure. 546 547 Collective on TS 548 549 Input Parameters: 550 + ts - the TS context obtained from TSCreate() 551 - viewer - visualization context 552 553 Options Database Key: 554 . -ts_view - calls TSView() at end of TSStep() 555 556 Notes: 557 The available visualization contexts include 558 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 559 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 560 output where only the first processor opens 561 the file. All other processors send their 562 data to the first processor to print. 563 564 The user can open an alternative visualization context with 565 PetscViewerASCIIOpen() - output to a specified file. 566 567 Level: beginner 568 569 .keywords: TS, timestep, view 570 571 .seealso: PetscViewerASCIIOpen() 572 @*/ 573 PetscErrorCode TSView(TS ts,PetscViewer viewer) 574 { 575 PetscErrorCode ierr; 576 char *type; 577 PetscTruth iascii,isstring; 578 579 PetscFunctionBegin; 580 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 581 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm); 582 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 583 PetscCheckSameComm(ts,1,viewer,2); 584 585 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 586 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 587 if (iascii) { 588 ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 589 ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr); 590 if (type) { 591 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 592 } else { 593 ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 594 } 595 if (ts->ops->view) { 596 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 597 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 598 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 599 } 600 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 601 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%g\n",ts->max_time);CHKERRQ(ierr); 602 if (ts->problem_type == TS_NONLINEAR) { 603 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 604 } 605 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 606 } else if (isstring) { 607 ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr); 608 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 609 } 610 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 611 if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 612 if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 613 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 614 PetscFunctionReturn(0); 615 } 616 617 618 #undef __FUNCT__ 619 #define __FUNCT__ "TSSetApplicationContext" 620 /*@C 621 TSSetApplicationContext - Sets an optional user-defined context for 622 the timesteppers. 623 624 Collective on TS 625 626 Input Parameters: 627 + ts - the TS context obtained from TSCreate() 628 - usrP - optional user context 629 630 Level: intermediate 631 632 .keywords: TS, timestep, set, application, context 633 634 .seealso: TSGetApplicationContext() 635 @*/ 636 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 637 { 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 640 ts->user = usrP; 641 PetscFunctionReturn(0); 642 } 643 644 #undef __FUNCT__ 645 #define __FUNCT__ "TSGetApplicationContext" 646 /*@C 647 TSGetApplicationContext - Gets the user-defined context for the 648 timestepper. 649 650 Not Collective 651 652 Input Parameter: 653 . ts - the TS context obtained from TSCreate() 654 655 Output Parameter: 656 . usrP - user context 657 658 Level: intermediate 659 660 .keywords: TS, timestep, get, application, context 661 662 .seealso: TSSetApplicationContext() 663 @*/ 664 PetscErrorCode TSGetApplicationContext(TS ts,void **usrP) 665 { 666 PetscFunctionBegin; 667 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 668 *usrP = ts->user; 669 PetscFunctionReturn(0); 670 } 671 672 #undef __FUNCT__ 673 #define __FUNCT__ "TSGetTimeStepNumber" 674 /*@ 675 TSGetTimeStepNumber - Gets the current number of timesteps. 676 677 Not Collective 678 679 Input Parameter: 680 . ts - the TS context obtained from TSCreate() 681 682 Output Parameter: 683 . iter - number steps so far 684 685 Level: intermediate 686 687 .keywords: TS, timestep, get, iteration, number 688 @*/ 689 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 690 { 691 PetscFunctionBegin; 692 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 693 PetscValidIntPointer(iter,2); 694 *iter = ts->steps; 695 PetscFunctionReturn(0); 696 } 697 698 #undef __FUNCT__ 699 #define __FUNCT__ "TSSetInitialTimeStep" 700 /*@ 701 TSSetInitialTimeStep - Sets the initial timestep to be used, 702 as well as the initial time. 703 704 Collective on TS 705 706 Input Parameters: 707 + ts - the TS context obtained from TSCreate() 708 . initial_time - the initial time 709 - time_step - the size of the timestep 710 711 Level: intermediate 712 713 .seealso: TSSetTimeStep(), TSGetTimeStep() 714 715 .keywords: TS, set, initial, timestep 716 @*/ 717 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 718 { 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 721 ts->time_step = time_step; 722 ts->initial_time_step = time_step; 723 ts->ptime = initial_time; 724 PetscFunctionReturn(0); 725 } 726 727 #undef __FUNCT__ 728 #define __FUNCT__ "TSSetTimeStep" 729 /*@ 730 TSSetTimeStep - Allows one to reset the timestep at any time, 731 useful for simple pseudo-timestepping codes. 732 733 Collective on TS 734 735 Input Parameters: 736 + ts - the TS context obtained from TSCreate() 737 - time_step - the size of the timestep 738 739 Level: intermediate 740 741 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 742 743 .keywords: TS, set, timestep 744 @*/ 745 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 746 { 747 PetscFunctionBegin; 748 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 749 ts->time_step = time_step; 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "TSGetTimeStep" 755 /*@ 756 TSGetTimeStep - Gets the current timestep size. 757 758 Not Collective 759 760 Input Parameter: 761 . ts - the TS context obtained from TSCreate() 762 763 Output Parameter: 764 . dt - the current timestep size 765 766 Level: intermediate 767 768 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 769 770 .keywords: TS, get, timestep 771 @*/ 772 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 773 { 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 776 PetscValidDoublePointer(dt,2); 777 *dt = ts->time_step; 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "TSGetSolution" 783 /*@C 784 TSGetSolution - Returns the solution at the present timestep. It 785 is valid to call this routine inside the function that you are evaluating 786 in order to move to the new timestep. This vector not changed until 787 the solution at the next timestep has been calculated. 788 789 Not Collective, but Vec returned is parallel if TS is parallel 790 791 Input Parameter: 792 . ts - the TS context obtained from TSCreate() 793 794 Output Parameter: 795 . v - the vector containing the solution 796 797 Level: intermediate 798 799 .seealso: TSGetTimeStep() 800 801 .keywords: TS, timestep, get, solution 802 @*/ 803 PetscErrorCode TSGetSolution(TS ts,Vec *v) 804 { 805 PetscFunctionBegin; 806 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 807 PetscValidPointer(v,2); 808 *v = ts->vec_sol_always; 809 PetscFunctionReturn(0); 810 } 811 812 /* ----- Routines to initialize and destroy a timestepper ---- */ 813 #undef __FUNCT__ 814 #define __FUNCT__ "TSSetProblemType" 815 /*@C 816 TSSetProblemType - Sets the type of problem to be solved. 817 818 Not collective 819 820 Input Parameters: 821 + ts - The TS 822 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 823 .vb 824 U_t = A U 825 U_t = A(t) U 826 U_t = F(t,U) 827 .ve 828 829 Level: beginner 830 831 .keywords: TS, problem type 832 .seealso: TSSetUp(), TSProblemType, TS 833 @*/ 834 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 835 { 836 PetscFunctionBegin; 837 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 838 ts->problem_type = type; 839 PetscFunctionReturn(0); 840 } 841 842 #undef __FUNCT__ 843 #define __FUNCT__ "TSGetProblemType" 844 /*@C 845 TSGetProblemType - Gets the type of problem to be solved. 846 847 Not collective 848 849 Input Parameter: 850 . ts - The TS 851 852 Output Parameter: 853 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 854 .vb 855 U_t = A U 856 U_t = A(t) U 857 U_t = F(t,U) 858 .ve 859 860 Level: beginner 861 862 .keywords: TS, problem type 863 .seealso: TSSetUp(), TSProblemType, TS 864 @*/ 865 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 866 { 867 PetscFunctionBegin; 868 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 869 PetscValidIntPointer(type,2); 870 *type = ts->problem_type; 871 PetscFunctionReturn(0); 872 } 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "TSSetUp" 876 /*@ 877 TSSetUp - Sets up the internal data structures for the later use 878 of a timestepper. 879 880 Collective on TS 881 882 Input Parameter: 883 . ts - the TS context obtained from TSCreate() 884 885 Notes: 886 For basic use of the TS solvers the user need not explicitly call 887 TSSetUp(), since these actions will automatically occur during 888 the call to TSStep(). However, if one wishes to control this 889 phase separately, TSSetUp() should be called after TSCreate() 890 and optional routines of the form TSSetXXX(), but before TSStep(). 891 892 Level: advanced 893 894 .keywords: TS, timestep, setup 895 896 .seealso: TSCreate(), TSStep(), TSDestroy() 897 @*/ 898 PetscErrorCode TSSetUp(TS ts) 899 { 900 PetscErrorCode ierr; 901 902 PetscFunctionBegin; 903 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 904 if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 905 if (!ts->type_name) { 906 ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 907 } 908 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 909 ts->setupcalled = 1; 910 PetscFunctionReturn(0); 911 } 912 913 #undef __FUNCT__ 914 #define __FUNCT__ "TSDestroy" 915 /*@C 916 TSDestroy - Destroys the timestepper context that was created 917 with TSCreate(). 918 919 Collective on TS 920 921 Input Parameter: 922 . ts - the TS context obtained from TSCreate() 923 924 Level: beginner 925 926 .keywords: TS, timestepper, destroy 927 928 .seealso: TSCreate(), TSSetUp(), TSSolve() 929 @*/ 930 PetscErrorCode TSDestroy(TS ts) 931 { 932 PetscErrorCode ierr; 933 PetscInt i; 934 935 PetscFunctionBegin; 936 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 937 if (--ts->refct > 0) PetscFunctionReturn(0); 938 939 /* if memory was published with AMS then destroy it */ 940 ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 941 942 if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 943 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 944 ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr); 945 for (i=0; i<ts->numbermonitors; i++) { 946 if (ts->mdestroy[i]) { 947 ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 948 } 949 } 950 PetscLogObjectDestroy((PetscObject)ts); 951 PetscHeaderDestroy((PetscObject)ts); 952 PetscFunctionReturn(0); 953 } 954 955 #undef __FUNCT__ 956 #define __FUNCT__ "TSGetSNES" 957 /*@C 958 TSGetSNES - Returns the SNES (nonlinear solver) associated with 959 a TS (timestepper) context. Valid only for nonlinear problems. 960 961 Not Collective, but SNES is parallel if TS is parallel 962 963 Input Parameter: 964 . ts - the TS context obtained from TSCreate() 965 966 Output Parameter: 967 . snes - the nonlinear solver context 968 969 Notes: 970 The user can then directly manipulate the SNES context to set various 971 options, etc. Likewise, the user can then extract and manipulate the 972 KSP, KSP, and PC contexts as well. 973 974 TSGetSNES() does not work for integrators that do not use SNES; in 975 this case TSGetSNES() returns PETSC_NULL in snes. 976 977 Level: beginner 978 979 .keywords: timestep, get, SNES 980 @*/ 981 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 982 { 983 PetscFunctionBegin; 984 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 985 PetscValidPointer(snes,2); 986 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 987 *snes = ts->snes; 988 PetscFunctionReturn(0); 989 } 990 991 #undef __FUNCT__ 992 #define __FUNCT__ "TSGetKSP" 993 /*@C 994 TSGetKSP - Returns the KSP (linear solver) associated with 995 a TS (timestepper) context. 996 997 Not Collective, but KSP is parallel if TS is parallel 998 999 Input Parameter: 1000 . ts - the TS context obtained from TSCreate() 1001 1002 Output Parameter: 1003 . ksp - the nonlinear solver context 1004 1005 Notes: 1006 The user can then directly manipulate the KSP context to set various 1007 options, etc. Likewise, the user can then extract and manipulate the 1008 KSP and PC contexts as well. 1009 1010 TSGetKSP() does not work for integrators that do not use KSP; 1011 in this case TSGetKSP() returns PETSC_NULL in ksp. 1012 1013 Level: beginner 1014 1015 .keywords: timestep, get, KSP 1016 @*/ 1017 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1018 { 1019 PetscFunctionBegin; 1020 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1021 PetscValidPointer(ksp,2); 1022 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1023 *ksp = ts->ksp; 1024 PetscFunctionReturn(0); 1025 } 1026 1027 /* ----------- Routines to set solver parameters ---------- */ 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "TSGetDuration" 1031 /*@ 1032 TSGetDuration - Gets the maximum number of timesteps to use and 1033 maximum time for iteration. 1034 1035 Collective on TS 1036 1037 Input Parameters: 1038 + ts - the TS context obtained from TSCreate() 1039 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1040 - maxtime - final time to iterate to, or PETSC_NULL 1041 1042 Level: intermediate 1043 1044 .keywords: TS, timestep, get, maximum, iterations, time 1045 @*/ 1046 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1047 { 1048 PetscFunctionBegin; 1049 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1050 if (maxsteps != PETSC_NULL) { 1051 PetscValidIntPointer(maxsteps,2); 1052 *maxsteps = ts->max_steps; 1053 } 1054 if (maxtime != PETSC_NULL) { 1055 PetscValidScalarPointer(maxtime,3); 1056 *maxtime = ts->max_time; 1057 } 1058 PetscFunctionReturn(0); 1059 } 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "TSSetDuration" 1063 /*@ 1064 TSSetDuration - Sets the maximum number of timesteps to use and 1065 maximum time for iteration. 1066 1067 Collective on TS 1068 1069 Input Parameters: 1070 + ts - the TS context obtained from TSCreate() 1071 . maxsteps - maximum number of iterations to use 1072 - maxtime - final time to iterate to 1073 1074 Options Database Keys: 1075 . -ts_max_steps <maxsteps> - Sets maxsteps 1076 . -ts_max_time <maxtime> - Sets maxtime 1077 1078 Notes: 1079 The default maximum number of iterations is 5000. Default time is 5.0 1080 1081 Level: intermediate 1082 1083 .keywords: TS, timestep, set, maximum, iterations 1084 @*/ 1085 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1086 { 1087 PetscFunctionBegin; 1088 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1089 ts->max_steps = maxsteps; 1090 ts->max_time = maxtime; 1091 PetscFunctionReturn(0); 1092 } 1093 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "TSSetSolution" 1096 /*@ 1097 TSSetSolution - Sets the initial solution vector 1098 for use by the TS routines. 1099 1100 Collective on TS and Vec 1101 1102 Input Parameters: 1103 + ts - the TS context obtained from TSCreate() 1104 - x - the solution vector 1105 1106 Level: beginner 1107 1108 .keywords: TS, timestep, set, solution, initial conditions 1109 @*/ 1110 PetscErrorCode TSSetSolution(TS ts,Vec x) 1111 { 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1114 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1115 ts->vec_sol = ts->vec_sol_always = x; 1116 PetscFunctionReturn(0); 1117 } 1118 1119 #undef __FUNCT__ 1120 #define __FUNCT__ "TSSetRhsBC" 1121 /*@C 1122 TSSetRhsBC - Sets the function which applies boundary conditions 1123 to the Rhs of each system. 1124 1125 Collective on TS 1126 1127 Input Parameters: 1128 + ts - The TS context obtained from TSCreate() 1129 - func - The function 1130 1131 Calling sequence of func: 1132 . func (TS ts, Vec rhs, void *ctx); 1133 1134 + rhs - The current rhs vector 1135 - ctx - The user-context 1136 1137 Level: intermediate 1138 1139 .keywords: TS, Rhs, boundary conditions 1140 @*/ 1141 PetscErrorCode TSSetRhsBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *)) 1142 { 1143 PetscFunctionBegin; 1144 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1145 ts->ops->applyrhsbc = func; 1146 PetscFunctionReturn(0); 1147 } 1148 1149 #undef __FUNCT__ 1150 #define __FUNCT__ "TSDefaultRhsBC" 1151 /*@ 1152 TSDefaultRhsBC - The default boundary condition function which does nothing. 1153 1154 Collective on TS 1155 1156 Input Parameters: 1157 + ts - The TS context obtained from TSCreate() 1158 . rhs - The Rhs 1159 - ctx - The user-context 1160 1161 Level: developer 1162 1163 .keywords: TS, Rhs, boundary conditions 1164 @*/ 1165 PetscErrorCode TSDefaultRhsBC(TS ts, Vec rhs, void *ctx) 1166 { 1167 PetscFunctionBegin; 1168 PetscFunctionReturn(0); 1169 } 1170 1171 #undef __FUNCT__ 1172 #define __FUNCT__ "TSSetSystemMatrixBC" 1173 /*@C 1174 TSSetSystemMatrixBC - Sets the function which applies boundary conditions 1175 to the system matrix and preconditioner of each system. 1176 1177 Collective on TS 1178 1179 Input Parameters: 1180 + ts - The TS context obtained from TSCreate() 1181 - func - The function 1182 1183 Calling sequence of func: 1184 . func (TS ts, Mat A, Mat B, void *ctx); 1185 1186 + A - The current system matrix 1187 . B - The current preconditioner 1188 - ctx - The user-context 1189 1190 Level: intermediate 1191 1192 .keywords: TS, System matrix, boundary conditions 1193 @*/ 1194 PetscErrorCode TSSetSystemMatrixBC(TS ts, PetscErrorCode (*func)(TS, Mat, Mat, void *)) 1195 { 1196 PetscFunctionBegin; 1197 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1198 ts->ops->applymatrixbc = func; 1199 PetscFunctionReturn(0); 1200 } 1201 1202 #undef __FUNCT__ 1203 #define __FUNCT__ "TSDefaultSystemMatrixBC" 1204 /*@ 1205 TSDefaultSystemMatrixBC - The default boundary condition function which 1206 does nothing. 1207 1208 Collective on TS 1209 1210 Input Parameters: 1211 + ts - The TS context obtained from TSCreate() 1212 . A - The system matrix 1213 . B - The preconditioner 1214 - ctx - The user-context 1215 1216 Level: developer 1217 1218 .keywords: TS, System matrix, boundary conditions 1219 @*/ 1220 PetscErrorCode TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx) 1221 { 1222 PetscFunctionBegin; 1223 PetscFunctionReturn(0); 1224 } 1225 1226 #undef __FUNCT__ 1227 #define __FUNCT__ "TSSetSolutionBC" 1228 /*@C 1229 TSSetSolutionBC - Sets the function which applies boundary conditions 1230 to the solution of each system. This is necessary in nonlinear systems 1231 which time dependent boundary conditions. 1232 1233 Collective on TS 1234 1235 Input Parameters: 1236 + ts - The TS context obtained from TSCreate() 1237 - func - The function 1238 1239 Calling sequence of func: 1240 . func (TS ts, Vec rsol, void *ctx); 1241 1242 + sol - The current solution vector 1243 - ctx - The user-context 1244 1245 Level: intermediate 1246 1247 .keywords: TS, solution, boundary conditions 1248 @*/ 1249 PetscErrorCode TSSetSolutionBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *)) 1250 { 1251 PetscFunctionBegin; 1252 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1253 ts->ops->applysolbc = func; 1254 PetscFunctionReturn(0); 1255 } 1256 1257 #undef __FUNCT__ 1258 #define __FUNCT__ "TSDefaultSolutionBC" 1259 /*@ 1260 TSDefaultSolutionBC - The default boundary condition function which 1261 does nothing. 1262 1263 Collective on TS 1264 1265 Input Parameters: 1266 + ts - The TS context obtained from TSCreate() 1267 . sol - The solution 1268 - ctx - The user-context 1269 1270 Level: developer 1271 1272 .keywords: TS, solution, boundary conditions 1273 @*/ 1274 PetscErrorCode TSDefaultSolutionBC(TS ts, Vec sol, void *ctx) 1275 { 1276 PetscFunctionBegin; 1277 PetscFunctionReturn(0); 1278 } 1279 1280 #undef __FUNCT__ 1281 #define __FUNCT__ "TSSetPreStep" 1282 /*@C 1283 TSSetPreStep - Sets the general-purpose function 1284 called once at the beginning of time stepping. 1285 1286 Collective on TS 1287 1288 Input Parameters: 1289 + ts - The TS context obtained from TSCreate() 1290 - func - The function 1291 1292 Calling sequence of func: 1293 . func (TS ts); 1294 1295 Level: intermediate 1296 1297 .keywords: TS, timestep 1298 @*/ 1299 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1300 { 1301 PetscFunctionBegin; 1302 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1303 ts->ops->prestep = func; 1304 PetscFunctionReturn(0); 1305 } 1306 1307 #undef __FUNCT__ 1308 #define __FUNCT__ "TSDefaultPreStep" 1309 /*@ 1310 TSDefaultPreStep - The default pre-stepping function which does nothing. 1311 1312 Collective on TS 1313 1314 Input Parameters: 1315 . ts - The TS context obtained from TSCreate() 1316 1317 Level: developer 1318 1319 .keywords: TS, timestep 1320 @*/ 1321 PetscErrorCode TSDefaultPreStep(TS ts) 1322 { 1323 PetscFunctionBegin; 1324 PetscFunctionReturn(0); 1325 } 1326 1327 #undef __FUNCT__ 1328 #define __FUNCT__ "TSSetUpdate" 1329 /*@C 1330 TSSetUpdate - Sets the general-purpose update function called 1331 at the beginning of every time step. This function can change 1332 the time step. 1333 1334 Collective on TS 1335 1336 Input Parameters: 1337 + ts - The TS context obtained from TSCreate() 1338 - func - The function 1339 1340 Calling sequence of func: 1341 . func (TS ts, double t, double *dt); 1342 1343 + t - The current time 1344 - dt - The current time step 1345 1346 Level: intermediate 1347 1348 .keywords: TS, update, timestep 1349 @*/ 1350 PetscErrorCode TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *)) 1351 { 1352 PetscFunctionBegin; 1353 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1354 ts->ops->update = func; 1355 PetscFunctionReturn(0); 1356 } 1357 1358 #undef __FUNCT__ 1359 #define __FUNCT__ "TSDefaultUpdate" 1360 /*@ 1361 TSDefaultUpdate - The default update function which does nothing. 1362 1363 Collective on TS 1364 1365 Input Parameters: 1366 + ts - The TS context obtained from TSCreate() 1367 - t - The current time 1368 1369 Output Parameters: 1370 . dt - The current time step 1371 1372 Level: developer 1373 1374 .keywords: TS, update, timestep 1375 @*/ 1376 PetscErrorCode TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt) 1377 { 1378 PetscFunctionBegin; 1379 PetscFunctionReturn(0); 1380 } 1381 1382 #undef __FUNCT__ 1383 #define __FUNCT__ "TSSetPostStep" 1384 /*@C 1385 TSSetPostStep - Sets the general-purpose function 1386 called once at the end of time stepping. 1387 1388 Collective on TS 1389 1390 Input Parameters: 1391 + ts - The TS context obtained from TSCreate() 1392 - func - The function 1393 1394 Calling sequence of func: 1395 . func (TS ts); 1396 1397 Level: intermediate 1398 1399 .keywords: TS, timestep 1400 @*/ 1401 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1402 { 1403 PetscFunctionBegin; 1404 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1405 ts->ops->poststep = func; 1406 PetscFunctionReturn(0); 1407 } 1408 1409 #undef __FUNCT__ 1410 #define __FUNCT__ "TSDefaultPostStep" 1411 /*@ 1412 TSDefaultPostStep - The default post-stepping function which does nothing. 1413 1414 Collective on TS 1415 1416 Input Parameters: 1417 . ts - The TS context obtained from TSCreate() 1418 1419 Level: developer 1420 1421 .keywords: TS, timestep 1422 @*/ 1423 PetscErrorCode TSDefaultPostStep(TS ts) 1424 { 1425 PetscFunctionBegin; 1426 PetscFunctionReturn(0); 1427 } 1428 1429 /* ------------ Routines to set performance monitoring options ----------- */ 1430 1431 #undef __FUNCT__ 1432 #define __FUNCT__ "TSSetMonitor" 1433 /*@C 1434 TSSetMonitor - Sets an ADDITIONAL function that is to be used at every 1435 timestep to display the iteration's progress. 1436 1437 Collective on TS 1438 1439 Input Parameters: 1440 + ts - the TS context obtained from TSCreate() 1441 . func - monitoring routine 1442 . mctx - [optional] user-defined context for private data for the 1443 monitor routine (use PETSC_NULL if no context is desired) 1444 - monitordestroy - [optional] routine that frees monitor context 1445 (may be PETSC_NULL) 1446 1447 Calling sequence of func: 1448 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1449 1450 + ts - the TS context 1451 . steps - iteration number 1452 . time - current time 1453 . x - current iterate 1454 - mctx - [optional] monitoring context 1455 1456 Notes: 1457 This routine adds an additional monitor to the list of monitors that 1458 already has been loaded. 1459 1460 Level: intermediate 1461 1462 .keywords: TS, timestep, set, monitor 1463 1464 .seealso: TSDefaultMonitor(), TSClearMonitor() 1465 @*/ 1466 PetscErrorCode TSSetMonitor(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1467 { 1468 PetscFunctionBegin; 1469 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1470 if (ts->numbermonitors >= MAXTSMONITORS) { 1471 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1472 } 1473 ts->monitor[ts->numbermonitors] = monitor; 1474 ts->mdestroy[ts->numbermonitors] = mdestroy; 1475 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1476 PetscFunctionReturn(0); 1477 } 1478 1479 #undef __FUNCT__ 1480 #define __FUNCT__ "TSClearMonitor" 1481 /*@C 1482 TSClearMonitor - Clears all the monitors that have been set on a time-step object. 1483 1484 Collective on TS 1485 1486 Input Parameters: 1487 . ts - the TS context obtained from TSCreate() 1488 1489 Notes: 1490 There is no way to remove a single, specific monitor. 1491 1492 Level: intermediate 1493 1494 .keywords: TS, timestep, set, monitor 1495 1496 .seealso: TSDefaultMonitor(), TSSetMonitor() 1497 @*/ 1498 PetscErrorCode TSClearMonitor(TS ts) 1499 { 1500 PetscFunctionBegin; 1501 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1502 ts->numbermonitors = 0; 1503 PetscFunctionReturn(0); 1504 } 1505 1506 #undef __FUNCT__ 1507 #define __FUNCT__ "TSDefaultMonitor" 1508 PetscErrorCode TSDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1509 { 1510 PetscErrorCode ierr; 1511 1512 PetscFunctionBegin; 1513 ierr = PetscPrintf(ts->comm,"timestep %D dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1514 PetscFunctionReturn(0); 1515 } 1516 1517 #undef __FUNCT__ 1518 #define __FUNCT__ "TSStep" 1519 /*@ 1520 TSStep - Steps the requested number of timesteps. 1521 1522 Collective on TS 1523 1524 Input Parameter: 1525 . ts - the TS context obtained from TSCreate() 1526 1527 Output Parameters: 1528 + steps - number of iterations until termination 1529 - ptime - time until termination 1530 1531 Level: beginner 1532 1533 .keywords: TS, timestep, solve 1534 1535 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1536 @*/ 1537 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1538 { 1539 PetscErrorCode ierr; 1540 1541 PetscFunctionBegin; 1542 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1543 if (!ts->setupcalled) { 1544 ierr = TSSetUp(ts);CHKERRQ(ierr); 1545 } 1546 1547 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1548 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1549 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1550 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1551 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1552 1553 if (!PetscPreLoadingOn) { 1554 ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr); 1555 } 1556 PetscFunctionReturn(0); 1557 } 1558 1559 #undef __FUNCT__ 1560 #define __FUNCT__ "TSMonitor" 1561 /* 1562 Runs the user provided monitor routines, if they exists. 1563 */ 1564 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1565 { 1566 PetscErrorCode ierr; 1567 PetscInt i,n = ts->numbermonitors; 1568 1569 PetscFunctionBegin; 1570 for (i=0; i<n; i++) { 1571 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1572 } 1573 PetscFunctionReturn(0); 1574 } 1575 1576 /* ------------------------------------------------------------------------*/ 1577 1578 #undef __FUNCT__ 1579 #define __FUNCT__ "TSLGMonitorCreate" 1580 /*@C 1581 TSLGMonitorCreate - Creates a line graph context for use with 1582 TS to monitor convergence of preconditioned residual norms. 1583 1584 Collective on TS 1585 1586 Input Parameters: 1587 + host - the X display to open, or null for the local machine 1588 . label - the title to put in the title bar 1589 . x, y - the screen coordinates of the upper left coordinate of the window 1590 - m, n - the screen width and height in pixels 1591 1592 Output Parameter: 1593 . draw - the drawing context 1594 1595 Options Database Key: 1596 . -ts_xmonitor - automatically sets line graph monitor 1597 1598 Notes: 1599 Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1600 1601 Level: intermediate 1602 1603 .keywords: TS, monitor, line graph, residual, seealso 1604 1605 .seealso: TSLGMonitorDestroy(), TSSetMonitor() 1606 1607 @*/ 1608 PetscErrorCode TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1609 { 1610 PetscDraw win; 1611 PetscErrorCode ierr; 1612 1613 PetscFunctionBegin; 1614 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1615 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1616 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1617 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1618 1619 PetscLogObjectParent(*draw,win); 1620 PetscFunctionReturn(0); 1621 } 1622 1623 #undef __FUNCT__ 1624 #define __FUNCT__ "TSLGMonitor" 1625 PetscErrorCode TSLGMonitor(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1626 { 1627 PetscDrawLG lg = (PetscDrawLG) monctx; 1628 PetscReal x,y = ptime; 1629 PetscErrorCode ierr; 1630 1631 PetscFunctionBegin; 1632 if (!monctx) { 1633 MPI_Comm comm; 1634 PetscViewer viewer; 1635 1636 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1637 viewer = PETSC_VIEWER_DRAW_(comm); 1638 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1639 } 1640 1641 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1642 x = (PetscReal)n; 1643 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1644 if (n < 20 || (n % 5)) { 1645 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1646 } 1647 PetscFunctionReturn(0); 1648 } 1649 1650 #undef __FUNCT__ 1651 #define __FUNCT__ "TSLGMonitorDestroy" 1652 /*@C 1653 TSLGMonitorDestroy - Destroys a line graph context that was created 1654 with TSLGMonitorCreate(). 1655 1656 Collective on PetscDrawLG 1657 1658 Input Parameter: 1659 . draw - the drawing context 1660 1661 Level: intermediate 1662 1663 .keywords: TS, monitor, line graph, destroy 1664 1665 .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor(); 1666 @*/ 1667 PetscErrorCode TSLGMonitorDestroy(PetscDrawLG drawlg) 1668 { 1669 PetscDraw draw; 1670 PetscErrorCode ierr; 1671 1672 PetscFunctionBegin; 1673 ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1674 ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1675 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1676 PetscFunctionReturn(0); 1677 } 1678 1679 #undef __FUNCT__ 1680 #define __FUNCT__ "TSGetTime" 1681 /*@ 1682 TSGetTime - Gets the current time. 1683 1684 Not Collective 1685 1686 Input Parameter: 1687 . ts - the TS context obtained from TSCreate() 1688 1689 Output Parameter: 1690 . t - the current time 1691 1692 Contributed by: Matthew Knepley 1693 1694 Level: beginner 1695 1696 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1697 1698 .keywords: TS, get, time 1699 @*/ 1700 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1701 { 1702 PetscFunctionBegin; 1703 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1704 PetscValidDoublePointer(t,2); 1705 *t = ts->ptime; 1706 PetscFunctionReturn(0); 1707 } 1708 1709 #undef __FUNCT__ 1710 #define __FUNCT__ "TSSetOptionsPrefix" 1711 /*@C 1712 TSSetOptionsPrefix - Sets the prefix used for searching for all 1713 TS options in the database. 1714 1715 Collective on TS 1716 1717 Input Parameter: 1718 + ts - The TS context 1719 - prefix - The prefix to prepend to all option names 1720 1721 Notes: 1722 A hyphen (-) must NOT be given at the beginning of the prefix name. 1723 The first character of all runtime options is AUTOMATICALLY the 1724 hyphen. 1725 1726 Contributed by: Matthew Knepley 1727 1728 Level: advanced 1729 1730 .keywords: TS, set, options, prefix, database 1731 1732 .seealso: TSSetFromOptions() 1733 1734 @*/ 1735 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1736 { 1737 PetscErrorCode ierr; 1738 1739 PetscFunctionBegin; 1740 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1741 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1742 switch(ts->problem_type) { 1743 case TS_NONLINEAR: 1744 ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1745 break; 1746 case TS_LINEAR: 1747 ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1748 break; 1749 } 1750 PetscFunctionReturn(0); 1751 } 1752 1753 1754 #undef __FUNCT__ 1755 #define __FUNCT__ "TSAppendOptionsPrefix" 1756 /*@C 1757 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1758 TS options in the database. 1759 1760 Collective on TS 1761 1762 Input Parameter: 1763 + ts - The TS context 1764 - prefix - The prefix to prepend to all option names 1765 1766 Notes: 1767 A hyphen (-) must NOT be given at the beginning of the prefix name. 1768 The first character of all runtime options is AUTOMATICALLY the 1769 hyphen. 1770 1771 Contributed by: Matthew Knepley 1772 1773 Level: advanced 1774 1775 .keywords: TS, append, options, prefix, database 1776 1777 .seealso: TSGetOptionsPrefix() 1778 1779 @*/ 1780 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1781 { 1782 PetscErrorCode ierr; 1783 1784 PetscFunctionBegin; 1785 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1786 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1787 switch(ts->problem_type) { 1788 case TS_NONLINEAR: 1789 ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1790 break; 1791 case TS_LINEAR: 1792 ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1793 break; 1794 } 1795 PetscFunctionReturn(0); 1796 } 1797 1798 #undef __FUNCT__ 1799 #define __FUNCT__ "TSGetOptionsPrefix" 1800 /*@C 1801 TSGetOptionsPrefix - Sets the prefix used for searching for all 1802 TS options in the database. 1803 1804 Not Collective 1805 1806 Input Parameter: 1807 . ts - The TS context 1808 1809 Output Parameter: 1810 . prefix - A pointer to the prefix string used 1811 1812 Contributed by: Matthew Knepley 1813 1814 Notes: On the fortran side, the user should pass in a string 'prifix' of 1815 sufficient length to hold the prefix. 1816 1817 Level: intermediate 1818 1819 .keywords: TS, get, options, prefix, database 1820 1821 .seealso: TSAppendOptionsPrefix() 1822 @*/ 1823 PetscErrorCode TSGetOptionsPrefix(TS ts,char *prefix[]) 1824 { 1825 PetscErrorCode ierr; 1826 1827 PetscFunctionBegin; 1828 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1829 PetscValidPointer(prefix,2); 1830 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1831 PetscFunctionReturn(0); 1832 } 1833 1834 #undef __FUNCT__ 1835 #define __FUNCT__ "TSGetRHSMatrix" 1836 /*@C 1837 TSGetRHSMatrix - Returns the matrix A at the present timestep. 1838 1839 Not Collective, but parallel objects are returned if TS is parallel 1840 1841 Input Parameter: 1842 . ts - The TS context obtained from TSCreate() 1843 1844 Output Parameters: 1845 + A - The matrix A, where U_t = A(t) U 1846 . M - The preconditioner matrix, usually the same as A 1847 - ctx - User-defined context for matrix evaluation routine 1848 1849 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1850 1851 Contributed by: Matthew Knepley 1852 1853 Level: intermediate 1854 1855 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 1856 1857 .keywords: TS, timestep, get, matrix 1858 1859 @*/ 1860 PetscErrorCode TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx) 1861 { 1862 PetscFunctionBegin; 1863 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1864 if (A) *A = ts->A; 1865 if (M) *M = ts->B; 1866 if (ctx) *ctx = ts->jacP; 1867 PetscFunctionReturn(0); 1868 } 1869 1870 #undef __FUNCT__ 1871 #define __FUNCT__ "TSGetRHSJacobian" 1872 /*@C 1873 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1874 1875 Not Collective, but parallel objects are returned if TS is parallel 1876 1877 Input Parameter: 1878 . ts - The TS context obtained from TSCreate() 1879 1880 Output Parameters: 1881 + J - The Jacobian J of F, where U_t = F(U,t) 1882 . M - The preconditioner matrix, usually the same as J 1883 - ctx - User-defined context for Jacobian evaluation routine 1884 1885 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1886 1887 Contributed by: Matthew Knepley 1888 1889 Level: intermediate 1890 1891 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber() 1892 1893 .keywords: TS, timestep, get, matrix, Jacobian 1894 @*/ 1895 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 1896 { 1897 PetscErrorCode ierr; 1898 1899 PetscFunctionBegin; 1900 ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr); 1901 PetscFunctionReturn(0); 1902 } 1903 1904 #undef __FUNCT__ 1905 #define __FUNCT__ "TSVecViewMonitor" 1906 /*@C 1907 TSVecViewMonitor - Monitors progress of the TS solvers by calling 1908 VecView() for the solution at each timestep 1909 1910 Collective on TS 1911 1912 Input Parameters: 1913 + ts - the TS context 1914 . step - current time-step 1915 . ptime - current time 1916 - dummy - either a viewer or PETSC_NULL 1917 1918 Level: intermediate 1919 1920 .keywords: TS, vector, monitor, view 1921 1922 .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView() 1923 @*/ 1924 PetscErrorCode TSVecViewMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 1925 { 1926 PetscErrorCode ierr; 1927 PetscViewer viewer = (PetscViewer) dummy; 1928 1929 PetscFunctionBegin; 1930 if (!viewer) { 1931 MPI_Comm comm; 1932 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1933 viewer = PETSC_VIEWER_DRAW_(comm); 1934 } 1935 ierr = VecView(x,viewer);CHKERRQ(ierr); 1936 PetscFunctionReturn(0); 1937 } 1938 1939 1940 1941