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