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((PetscObject)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__ "TSSetRhsBC" 1125 /*@C 1126 TSSetRhsBC - Sets the function which applies boundary conditions 1127 to the Rhs of each system. 1128 1129 Collective on TS 1130 1131 Input Parameters: 1132 + ts - The TS context obtained from TSCreate() 1133 - func - The function 1134 1135 Calling sequence of func: 1136 . func (TS ts, Vec rhs, void *ctx); 1137 1138 + rhs - The current rhs vector 1139 - ctx - The user-context 1140 1141 Level: intermediate 1142 1143 .keywords: TS, Rhs, boundary conditions 1144 @*/ 1145 PetscErrorCode TSSetRhsBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *)) 1146 { 1147 PetscFunctionBegin; 1148 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1149 ts->ops->applyrhsbc = func; 1150 PetscFunctionReturn(0); 1151 } 1152 1153 #undef __FUNCT__ 1154 #define __FUNCT__ "TSDefaultRhsBC" 1155 /*@ 1156 TSDefaultRhsBC - The default boundary condition function which does nothing. 1157 1158 Collective on TS 1159 1160 Input Parameters: 1161 + ts - The TS context obtained from TSCreate() 1162 . rhs - The Rhs 1163 - ctx - The user-context 1164 1165 Level: developer 1166 1167 .keywords: TS, Rhs, boundary conditions 1168 @*/ 1169 PetscErrorCode TSDefaultRhsBC(TS ts, Vec rhs, void *ctx) 1170 { 1171 PetscFunctionBegin; 1172 PetscFunctionReturn(0); 1173 } 1174 1175 #undef __FUNCT__ 1176 #define __FUNCT__ "TSSetSystemMatrixBC" 1177 /*@C 1178 TSSetSystemMatrixBC - Sets the function which applies boundary conditions 1179 to the system matrix and preconditioner of each system. 1180 1181 Collective on TS 1182 1183 Input Parameters: 1184 + ts - The TS context obtained from TSCreate() 1185 - func - The function 1186 1187 Calling sequence of func: 1188 . func (TS ts, Mat A, Mat B, void *ctx); 1189 1190 + A - The current system matrix 1191 . B - The current preconditioner 1192 - ctx - The user-context 1193 1194 Level: intermediate 1195 1196 .keywords: TS, System matrix, boundary conditions 1197 @*/ 1198 PetscErrorCode TSSetSystemMatrixBC(TS ts, PetscErrorCode (*func)(TS, Mat, Mat, void *)) 1199 { 1200 PetscFunctionBegin; 1201 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1202 ts->ops->applymatrixbc = func; 1203 PetscFunctionReturn(0); 1204 } 1205 1206 #undef __FUNCT__ 1207 #define __FUNCT__ "TSDefaultSystemMatrixBC" 1208 /*@ 1209 TSDefaultSystemMatrixBC - The default boundary condition function which 1210 does nothing. 1211 1212 Collective on TS 1213 1214 Input Parameters: 1215 + ts - The TS context obtained from TSCreate() 1216 . A - The system matrix 1217 . B - The preconditioner 1218 - ctx - The user-context 1219 1220 Level: developer 1221 1222 .keywords: TS, System matrix, boundary conditions 1223 @*/ 1224 PetscErrorCode TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx) 1225 { 1226 PetscFunctionBegin; 1227 PetscFunctionReturn(0); 1228 } 1229 1230 #undef __FUNCT__ 1231 #define __FUNCT__ "TSSetSolutionBC" 1232 /*@C 1233 TSSetSolutionBC - Sets the function which applies boundary conditions 1234 to the solution of each system. This is necessary in nonlinear systems 1235 which time dependent boundary conditions. 1236 1237 Collective on TS 1238 1239 Input Parameters: 1240 + ts - The TS context obtained from TSCreate() 1241 - func - The function 1242 1243 Calling sequence of func: 1244 . func (TS ts, Vec rsol, void *ctx); 1245 1246 + sol - The current solution vector 1247 - ctx - The user-context 1248 1249 Level: intermediate 1250 1251 .keywords: TS, solution, boundary conditions 1252 @*/ 1253 PetscErrorCode TSSetSolutionBC(TS ts, PetscErrorCode (*func)(TS, Vec, void *)) 1254 { 1255 PetscFunctionBegin; 1256 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1257 ts->ops->applysolbc = func; 1258 PetscFunctionReturn(0); 1259 } 1260 1261 #undef __FUNCT__ 1262 #define __FUNCT__ "TSDefaultSolutionBC" 1263 /*@ 1264 TSDefaultSolutionBC - The default boundary condition function which 1265 does nothing. 1266 1267 Collective on TS 1268 1269 Input Parameters: 1270 + ts - The TS context obtained from TSCreate() 1271 . sol - The solution 1272 - ctx - The user-context 1273 1274 Level: developer 1275 1276 .keywords: TS, solution, boundary conditions 1277 @*/ 1278 PetscErrorCode TSDefaultSolutionBC(TS ts, Vec sol, void *ctx) 1279 { 1280 PetscFunctionBegin; 1281 PetscFunctionReturn(0); 1282 } 1283 1284 #undef __FUNCT__ 1285 #define __FUNCT__ "TSSetPreStep" 1286 /*@C 1287 TSSetPreStep - Sets the general-purpose function 1288 called once at the beginning of time stepping. 1289 1290 Collective on TS 1291 1292 Input Parameters: 1293 + ts - The TS context obtained from TSCreate() 1294 - func - The function 1295 1296 Calling sequence of func: 1297 . func (TS ts); 1298 1299 Level: intermediate 1300 1301 .keywords: TS, timestep 1302 @*/ 1303 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1304 { 1305 PetscFunctionBegin; 1306 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1307 ts->ops->prestep = func; 1308 PetscFunctionReturn(0); 1309 } 1310 1311 #undef __FUNCT__ 1312 #define __FUNCT__ "TSDefaultPreStep" 1313 /*@ 1314 TSDefaultPreStep - The default pre-stepping function which does nothing. 1315 1316 Collective on TS 1317 1318 Input Parameters: 1319 . ts - The TS context obtained from TSCreate() 1320 1321 Level: developer 1322 1323 .keywords: TS, timestep 1324 @*/ 1325 PetscErrorCode TSDefaultPreStep(TS ts) 1326 { 1327 PetscFunctionBegin; 1328 PetscFunctionReturn(0); 1329 } 1330 1331 #undef __FUNCT__ 1332 #define __FUNCT__ "TSSetUpdate" 1333 /*@C 1334 TSSetUpdate - Sets the general-purpose update function called 1335 at the beginning of every time step. This function can change 1336 the time step. 1337 1338 Collective on TS 1339 1340 Input Parameters: 1341 + ts - The TS context obtained from TSCreate() 1342 - func - The function 1343 1344 Calling sequence of func: 1345 . func (TS ts, double t, double *dt); 1346 1347 + t - The current time 1348 - dt - The current time step 1349 1350 Level: intermediate 1351 1352 .keywords: TS, update, timestep 1353 @*/ 1354 PetscErrorCode TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *)) 1355 { 1356 PetscFunctionBegin; 1357 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1358 ts->ops->update = func; 1359 PetscFunctionReturn(0); 1360 } 1361 1362 #undef __FUNCT__ 1363 #define __FUNCT__ "TSDefaultUpdate" 1364 /*@ 1365 TSDefaultUpdate - The default update function which does nothing. 1366 1367 Collective on TS 1368 1369 Input Parameters: 1370 + ts - The TS context obtained from TSCreate() 1371 - t - The current time 1372 1373 Output Parameters: 1374 . dt - The current time step 1375 1376 Level: developer 1377 1378 .keywords: TS, update, timestep 1379 @*/ 1380 PetscErrorCode TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt) 1381 { 1382 PetscFunctionBegin; 1383 PetscFunctionReturn(0); 1384 } 1385 1386 #undef __FUNCT__ 1387 #define __FUNCT__ "TSSetPostStep" 1388 /*@C 1389 TSSetPostStep - Sets the general-purpose function 1390 called once at the end of time stepping. 1391 1392 Collective on TS 1393 1394 Input Parameters: 1395 + ts - The TS context obtained from TSCreate() 1396 - func - The function 1397 1398 Calling sequence of func: 1399 . func (TS ts); 1400 1401 Level: intermediate 1402 1403 .keywords: TS, timestep 1404 @*/ 1405 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1406 { 1407 PetscFunctionBegin; 1408 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1409 ts->ops->poststep = func; 1410 PetscFunctionReturn(0); 1411 } 1412 1413 #undef __FUNCT__ 1414 #define __FUNCT__ "TSDefaultPostStep" 1415 /*@ 1416 TSDefaultPostStep - The default post-stepping function which does nothing. 1417 1418 Collective on TS 1419 1420 Input Parameters: 1421 . ts - The TS context obtained from TSCreate() 1422 1423 Level: developer 1424 1425 .keywords: TS, timestep 1426 @*/ 1427 PetscErrorCode TSDefaultPostStep(TS ts) 1428 { 1429 PetscFunctionBegin; 1430 PetscFunctionReturn(0); 1431 } 1432 1433 /* ------------ Routines to set performance monitoring options ----------- */ 1434 1435 #undef __FUNCT__ 1436 #define __FUNCT__ "TSSetMonitor" 1437 /*@C 1438 TSSetMonitor - Sets an ADDITIONAL function that is to be used at every 1439 timestep to display the iteration's progress. 1440 1441 Collective on TS 1442 1443 Input Parameters: 1444 + ts - the TS context obtained from TSCreate() 1445 . func - monitoring routine 1446 . mctx - [optional] user-defined context for private data for the 1447 monitor routine (use PETSC_NULL if no context is desired) 1448 - monitordestroy - [optional] routine that frees monitor context 1449 (may be PETSC_NULL) 1450 1451 Calling sequence of func: 1452 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1453 1454 + ts - the TS context 1455 . steps - iteration number 1456 . time - current time 1457 . x - current iterate 1458 - mctx - [optional] monitoring context 1459 1460 Notes: 1461 This routine adds an additional monitor to the list of monitors that 1462 already has been loaded. 1463 1464 Level: intermediate 1465 1466 .keywords: TS, timestep, set, monitor 1467 1468 .seealso: TSDefaultMonitor(), TSClearMonitor() 1469 @*/ 1470 PetscErrorCode TSSetMonitor(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1471 { 1472 PetscFunctionBegin; 1473 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1474 if (ts->numbermonitors >= MAXTSMONITORS) { 1475 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1476 } 1477 ts->monitor[ts->numbermonitors] = monitor; 1478 ts->mdestroy[ts->numbermonitors] = mdestroy; 1479 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1480 PetscFunctionReturn(0); 1481 } 1482 1483 #undef __FUNCT__ 1484 #define __FUNCT__ "TSClearMonitor" 1485 /*@C 1486 TSClearMonitor - Clears all the monitors that have been set on a time-step object. 1487 1488 Collective on TS 1489 1490 Input Parameters: 1491 . ts - the TS context obtained from TSCreate() 1492 1493 Notes: 1494 There is no way to remove a single, specific monitor. 1495 1496 Level: intermediate 1497 1498 .keywords: TS, timestep, set, monitor 1499 1500 .seealso: TSDefaultMonitor(), TSSetMonitor() 1501 @*/ 1502 PetscErrorCode TSClearMonitor(TS ts) 1503 { 1504 PetscFunctionBegin; 1505 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1506 ts->numbermonitors = 0; 1507 PetscFunctionReturn(0); 1508 } 1509 1510 #undef __FUNCT__ 1511 #define __FUNCT__ "TSDefaultMonitor" 1512 PetscErrorCode TSDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1513 { 1514 PetscErrorCode ierr; 1515 1516 PetscFunctionBegin; 1517 ierr = PetscPrintf(ts->comm,"timestep %D dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1518 PetscFunctionReturn(0); 1519 } 1520 1521 #undef __FUNCT__ 1522 #define __FUNCT__ "TSStep" 1523 /*@ 1524 TSStep - Steps the requested number of timesteps. 1525 1526 Collective on TS 1527 1528 Input Parameter: 1529 . ts - the TS context obtained from TSCreate() 1530 1531 Output Parameters: 1532 + steps - number of iterations until termination 1533 - ptime - time until termination 1534 1535 Level: beginner 1536 1537 .keywords: TS, timestep, solve 1538 1539 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1540 @*/ 1541 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1542 { 1543 PetscErrorCode ierr; 1544 1545 PetscFunctionBegin; 1546 PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1547 if (!ts->setupcalled) { 1548 ierr = TSSetUp(ts);CHKERRQ(ierr); 1549 } 1550 1551 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1552 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1553 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1554 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1555 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1556 1557 if (!PetscPreLoadingOn) { 1558 ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr); 1559 } 1560 PetscFunctionReturn(0); 1561 } 1562 1563 #undef __FUNCT__ 1564 #define __FUNCT__ "TSMonitor" 1565 /* 1566 Runs the user provided monitor routines, if they exists. 1567 */ 1568 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1569 { 1570 PetscErrorCode ierr; 1571 PetscInt i,n = ts->numbermonitors; 1572 1573 PetscFunctionBegin; 1574 for (i=0; i<n; i++) { 1575 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1576 } 1577 PetscFunctionReturn(0); 1578 } 1579 1580 /* ------------------------------------------------------------------------*/ 1581 1582 #undef __FUNCT__ 1583 #define __FUNCT__ "TSLGMonitorCreate" 1584 /*@C 1585 TSLGMonitorCreate - Creates a line graph context for use with 1586 TS to monitor convergence of preconditioned residual norms. 1587 1588 Collective on TS 1589 1590 Input Parameters: 1591 + host - the X display to open, or null for the local machine 1592 . label - the title to put in the title bar 1593 . x, y - the screen coordinates of the upper left coordinate of the window 1594 - m, n - the screen width and height in pixels 1595 1596 Output Parameter: 1597 . draw - the drawing context 1598 1599 Options Database Key: 1600 . -ts_xmonitor - automatically sets line graph monitor 1601 1602 Notes: 1603 Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1604 1605 Level: intermediate 1606 1607 .keywords: TS, monitor, line graph, residual, seealso 1608 1609 .seealso: TSLGMonitorDestroy(), TSSetMonitor() 1610 1611 @*/ 1612 PetscErrorCode TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1613 { 1614 PetscDraw win; 1615 PetscErrorCode ierr; 1616 1617 PetscFunctionBegin; 1618 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1619 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1620 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1621 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1622 1623 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1624 PetscFunctionReturn(0); 1625 } 1626 1627 #undef __FUNCT__ 1628 #define __FUNCT__ "TSLGMonitor" 1629 PetscErrorCode TSLGMonitor(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1630 { 1631 PetscDrawLG lg = (PetscDrawLG) monctx; 1632 PetscReal x,y = ptime; 1633 PetscErrorCode ierr; 1634 1635 PetscFunctionBegin; 1636 if (!monctx) { 1637 MPI_Comm comm; 1638 PetscViewer viewer; 1639 1640 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1641 viewer = PETSC_VIEWER_DRAW_(comm); 1642 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1643 } 1644 1645 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1646 x = (PetscReal)n; 1647 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1648 if (n < 20 || (n % 5)) { 1649 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1650 } 1651 PetscFunctionReturn(0); 1652 } 1653 1654 #undef __FUNCT__ 1655 #define __FUNCT__ "TSLGMonitorDestroy" 1656 /*@C 1657 TSLGMonitorDestroy - Destroys a line graph context that was created 1658 with TSLGMonitorCreate(). 1659 1660 Collective on PetscDrawLG 1661 1662 Input Parameter: 1663 . draw - the drawing context 1664 1665 Level: intermediate 1666 1667 .keywords: TS, monitor, line graph, destroy 1668 1669 .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor(); 1670 @*/ 1671 PetscErrorCode TSLGMonitorDestroy(PetscDrawLG drawlg) 1672 { 1673 PetscDraw draw; 1674 PetscErrorCode ierr; 1675 1676 PetscFunctionBegin; 1677 ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1678 ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1679 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1680 PetscFunctionReturn(0); 1681 } 1682 1683 #undef __FUNCT__ 1684 #define __FUNCT__ "TSGetTime" 1685 /*@ 1686 TSGetTime - Gets the current time. 1687 1688 Not Collective 1689 1690 Input Parameter: 1691 . ts - the TS context obtained from TSCreate() 1692 1693 Output Parameter: 1694 . t - the current time 1695 1696 Contributed by: Matthew Knepley 1697 1698 Level: beginner 1699 1700 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1701 1702 .keywords: TS, get, time 1703 @*/ 1704 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1705 { 1706 PetscFunctionBegin; 1707 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1708 PetscValidDoublePointer(t,2); 1709 *t = ts->ptime; 1710 PetscFunctionReturn(0); 1711 } 1712 1713 #undef __FUNCT__ 1714 #define __FUNCT__ "TSSetOptionsPrefix" 1715 /*@C 1716 TSSetOptionsPrefix - Sets the prefix used for searching for all 1717 TS options in the database. 1718 1719 Collective on TS 1720 1721 Input Parameter: 1722 + ts - The TS context 1723 - prefix - The prefix to prepend to all option names 1724 1725 Notes: 1726 A hyphen (-) must NOT be given at the beginning of the prefix name. 1727 The first character of all runtime options is AUTOMATICALLY the 1728 hyphen. 1729 1730 Contributed by: Matthew Knepley 1731 1732 Level: advanced 1733 1734 .keywords: TS, set, options, prefix, database 1735 1736 .seealso: TSSetFromOptions() 1737 1738 @*/ 1739 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1740 { 1741 PetscErrorCode ierr; 1742 1743 PetscFunctionBegin; 1744 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1745 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1746 switch(ts->problem_type) { 1747 case TS_NONLINEAR: 1748 ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1749 break; 1750 case TS_LINEAR: 1751 ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1752 break; 1753 } 1754 PetscFunctionReturn(0); 1755 } 1756 1757 1758 #undef __FUNCT__ 1759 #define __FUNCT__ "TSAppendOptionsPrefix" 1760 /*@C 1761 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1762 TS options in the database. 1763 1764 Collective on TS 1765 1766 Input Parameter: 1767 + ts - The TS context 1768 - prefix - The prefix to prepend to all option names 1769 1770 Notes: 1771 A hyphen (-) must NOT be given at the beginning of the prefix name. 1772 The first character of all runtime options is AUTOMATICALLY the 1773 hyphen. 1774 1775 Contributed by: Matthew Knepley 1776 1777 Level: advanced 1778 1779 .keywords: TS, append, options, prefix, database 1780 1781 .seealso: TSGetOptionsPrefix() 1782 1783 @*/ 1784 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1785 { 1786 PetscErrorCode ierr; 1787 1788 PetscFunctionBegin; 1789 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1790 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1791 switch(ts->problem_type) { 1792 case TS_NONLINEAR: 1793 ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1794 break; 1795 case TS_LINEAR: 1796 ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1797 break; 1798 } 1799 PetscFunctionReturn(0); 1800 } 1801 1802 #undef __FUNCT__ 1803 #define __FUNCT__ "TSGetOptionsPrefix" 1804 /*@C 1805 TSGetOptionsPrefix - Sets the prefix used for searching for all 1806 TS options in the database. 1807 1808 Not Collective 1809 1810 Input Parameter: 1811 . ts - The TS context 1812 1813 Output Parameter: 1814 . prefix - A pointer to the prefix string used 1815 1816 Contributed by: Matthew Knepley 1817 1818 Notes: On the fortran side, the user should pass in a string 'prifix' of 1819 sufficient length to hold the prefix. 1820 1821 Level: intermediate 1822 1823 .keywords: TS, get, options, prefix, database 1824 1825 .seealso: TSAppendOptionsPrefix() 1826 @*/ 1827 PetscErrorCode TSGetOptionsPrefix(TS ts,char *prefix[]) 1828 { 1829 PetscErrorCode ierr; 1830 1831 PetscFunctionBegin; 1832 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1833 PetscValidPointer(prefix,2); 1834 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1835 PetscFunctionReturn(0); 1836 } 1837 1838 #undef __FUNCT__ 1839 #define __FUNCT__ "TSGetRHSMatrix" 1840 /*@C 1841 TSGetRHSMatrix - Returns the matrix A at the present timestep. 1842 1843 Not Collective, but parallel objects are returned if TS is parallel 1844 1845 Input Parameter: 1846 . ts - The TS context obtained from TSCreate() 1847 1848 Output Parameters: 1849 + A - The matrix A, where U_t = A(t) U 1850 . M - The preconditioner matrix, usually the same as A 1851 - ctx - User-defined context for matrix evaluation routine 1852 1853 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1854 1855 Contributed by: Matthew Knepley 1856 1857 Level: intermediate 1858 1859 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 1860 1861 .keywords: TS, timestep, get, matrix 1862 1863 @*/ 1864 PetscErrorCode TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx) 1865 { 1866 PetscFunctionBegin; 1867 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1868 if (A) *A = ts->A; 1869 if (M) *M = ts->B; 1870 if (ctx) *ctx = ts->jacP; 1871 PetscFunctionReturn(0); 1872 } 1873 1874 #undef __FUNCT__ 1875 #define __FUNCT__ "TSGetRHSJacobian" 1876 /*@C 1877 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1878 1879 Not Collective, but parallel objects are returned if TS is parallel 1880 1881 Input Parameter: 1882 . ts - The TS context obtained from TSCreate() 1883 1884 Output Parameters: 1885 + J - The Jacobian J of F, where U_t = F(U,t) 1886 . M - The preconditioner matrix, usually the same as J 1887 - ctx - User-defined context for Jacobian evaluation routine 1888 1889 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1890 1891 Contributed by: Matthew Knepley 1892 1893 Level: intermediate 1894 1895 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber() 1896 1897 .keywords: TS, timestep, get, matrix, Jacobian 1898 @*/ 1899 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 1900 { 1901 PetscErrorCode ierr; 1902 1903 PetscFunctionBegin; 1904 ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr); 1905 PetscFunctionReturn(0); 1906 } 1907 1908 #undef __FUNCT__ 1909 #define __FUNCT__ "TSVecViewMonitor" 1910 /*@C 1911 TSVecViewMonitor - Monitors progress of the TS solvers by calling 1912 VecView() for the solution at each timestep 1913 1914 Collective on TS 1915 1916 Input Parameters: 1917 + ts - the TS context 1918 . step - current time-step 1919 . ptime - current time 1920 - dummy - either a viewer or PETSC_NULL 1921 1922 Level: intermediate 1923 1924 .keywords: TS, vector, monitor, view 1925 1926 .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView() 1927 @*/ 1928 PetscErrorCode TSVecViewMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 1929 { 1930 PetscErrorCode ierr; 1931 PetscViewer viewer = (PetscViewer) dummy; 1932 1933 PetscFunctionBegin; 1934 if (!viewer) { 1935 MPI_Comm comm; 1936 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1937 viewer = PETSC_VIEWER_DRAW_(comm); 1938 } 1939 ierr = VecView(x,viewer);CHKERRQ(ierr); 1940 PetscFunctionReturn(0); 1941 } 1942 1943 1944 1945