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