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