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