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