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