1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 6 typedef struct { 7 Vec update; /* work vector where new solution is formed */ 8 Vec func; /* work vector where F(t[i],u[i]) is stored */ 9 Vec xdot; /* work vector for time derivative of state */ 10 11 /* information used for Pseudo-timestepping */ 12 13 PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */ 14 void *dtctx; 15 PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */ 16 void *verifyctx; 17 18 PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_initial; /* initial time-step */ 22 PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 23 PetscReal dt_max; /* maximum time step */ 24 PetscBool increment_dt_from_initial_dt; 25 PetscReal fatol, frtol; 26 } TS_Pseudo; 27 28 /* ------------------------------------------------------------------------------*/ 29 30 /*@C 31 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 32 pseudo-timestepping process. 33 34 Collective on ts 35 36 Input Parameter: 37 . ts - timestep context 38 39 Output Parameter: 40 . dt - newly computed timestep 41 42 Level: developer 43 44 Note: 45 The routine to be called here to compute the timestep should be 46 set by calling `TSPseudoSetTimeStep()`. 47 48 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 49 @*/ 50 PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 51 { 52 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53 54 PetscFunctionBegin; 55 PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 56 PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx)); 57 PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 58 PetscFunctionReturn(0); 59 } 60 61 /* ------------------------------------------------------------------------------*/ 62 /*@C 63 TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 64 65 Collective on ts 66 67 Input Parameters: 68 + ts - the timestep context 69 . dtctx - unused timestep context 70 - update - latest solution vector 71 72 Output Parameters: 73 + newdt - the timestep to use for the next step 74 - flag - flag indicating whether the last time step was acceptable 75 76 Level: advanced 77 78 Note: 79 This routine always returns a flag of 1, indicating an acceptable 80 timestep. 81 82 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 83 @*/ 84 PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 85 { 86 PetscFunctionBegin; 87 *flag = PETSC_TRUE; 88 PetscFunctionReturn(0); 89 } 90 91 /*@ 92 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 93 94 Collective on ts 95 96 Input Parameters: 97 + ts - timestep context 98 - update - latest solution vector 99 100 Output Parameters: 101 + dt - newly computed timestep (if it had to shrink) 102 - flag - indicates if current timestep was ok 103 104 Level: advanced 105 106 Notes: 107 The routine to be called here to compute the timestep should be 108 set by calling `TSPseudoSetVerifyTimeStep()`. 109 110 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 111 @*/ 112 PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 113 { 114 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 115 116 PetscFunctionBegin; 117 *flag = PETSC_TRUE; 118 if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag)); 119 PetscFunctionReturn(0); 120 } 121 122 /* --------------------------------------------------------------------------------*/ 123 124 static PetscErrorCode TSStep_Pseudo(TS ts) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 127 PetscInt nits, lits, reject; 128 PetscBool stepok; 129 PetscReal next_time_step = ts->time_step; 130 131 PetscFunctionBegin; 132 if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 133 PetscCall(VecCopy(ts->vec_sol, pseudo->update)); 134 PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step)); 135 for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) { 136 ts->time_step = next_time_step; 137 PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 138 PetscCall(SNESSolve(ts->snes, NULL, pseudo->update)); 139 PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 140 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 141 ts->snes_its += nits; 142 ts->ksp_its += lits; 143 PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &(pseudo->update))); 144 PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok)); 145 if (!stepok) { 146 next_time_step = ts->time_step; 147 continue; 148 } 149 pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 150 PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok)); 151 if (stepok) break; 152 } 153 if (reject >= ts->max_reject) { 154 ts->reason = TS_DIVERGED_STEP_REJECTED; 155 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject)); 156 PetscFunctionReturn(0); 157 } 158 159 PetscCall(VecCopy(pseudo->update, ts->vec_sol)); 160 ts->ptime += ts->time_step; 161 ts->time_step = next_time_step; 162 163 if (pseudo->fnorm < 0) { 164 PetscCall(VecZeroEntries(pseudo->xdot)); 165 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 166 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 167 } 168 if (pseudo->fnorm < pseudo->fatol) { 169 ts->reason = TS_CONVERGED_PSEUDO_FATOL; 170 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 171 PetscFunctionReturn(0); 172 } 173 if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) { 174 ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 175 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g / fnorm_initial %g < frtol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->fnorm_initial, (double)pseudo->fatol)); 176 PetscFunctionReturn(0); 177 } 178 PetscFunctionReturn(0); 179 } 180 181 /*------------------------------------------------------------*/ 182 static PetscErrorCode TSReset_Pseudo(TS ts) 183 { 184 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 185 186 PetscFunctionBegin; 187 PetscCall(VecDestroy(&pseudo->update)); 188 PetscCall(VecDestroy(&pseudo->func)); 189 PetscCall(VecDestroy(&pseudo->xdot)); 190 PetscFunctionReturn(0); 191 } 192 193 static PetscErrorCode TSDestroy_Pseudo(TS ts) 194 { 195 PetscFunctionBegin; 196 PetscCall(TSReset_Pseudo(ts)); 197 PetscCall(PetscFree(ts->data)); 198 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 199 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 200 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 201 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 202 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 203 PetscFunctionReturn(0); 204 } 205 206 /*------------------------------------------------------------*/ 207 208 /* 209 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 210 */ 211 static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 212 { 213 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 214 const PetscScalar mdt = 1.0 / ts->time_step, *xnp1, *xn; 215 PetscScalar *xdot; 216 PetscInt i, n; 217 218 PetscFunctionBegin; 219 *Xdot = NULL; 220 PetscCall(VecGetArrayRead(ts->vec_sol, &xn)); 221 PetscCall(VecGetArrayRead(X, &xnp1)); 222 PetscCall(VecGetArray(pseudo->xdot, &xdot)); 223 PetscCall(VecGetLocalSize(X, &n)); 224 for (i = 0; i < n; i++) xdot[i] = mdt * (xnp1[i] - xn[i]); 225 PetscCall(VecRestoreArrayRead(ts->vec_sol, &xn)); 226 PetscCall(VecRestoreArrayRead(X, &xnp1)); 227 PetscCall(VecRestoreArray(pseudo->xdot, &xdot)); 228 *Xdot = pseudo->xdot; 229 PetscFunctionReturn(0); 230 } 231 232 /* 233 The transient residual is 234 235 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 236 237 or for ODE, 238 239 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 240 241 This is the function that must be evaluated for transient simulation and for 242 finite difference Jacobians. On the first Newton step, this algorithm uses 243 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 244 residual is actually the steady state residual. Pseudotransient 245 continuation as described in the literature is a linearly implicit 246 algorithm, it only takes this one Newton step with the steady state 247 residual, and then advances to the next time step. 248 */ 249 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 250 { 251 Vec Xdot; 252 253 PetscFunctionBegin; 254 PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 255 PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, Xdot, Y, PETSC_FALSE)); 256 PetscFunctionReturn(0); 257 } 258 259 /* 260 This constructs the Jacobian needed for SNES. For DAE, this is 261 262 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 263 264 and for ODE: 265 266 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 267 */ 268 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 269 { 270 Vec Xdot; 271 272 PetscFunctionBegin; 273 PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 274 PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode TSSetUp_Pseudo(TS ts) 279 { 280 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 281 282 PetscFunctionBegin; 283 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update)); 284 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 285 PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 286 PetscFunctionReturn(0); 287 } 288 /*------------------------------------------------------------*/ 289 290 static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 291 { 292 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 293 PetscViewer viewer = (PetscViewer)dummy; 294 295 PetscFunctionBegin; 296 if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 297 PetscCall(VecZeroEntries(pseudo->xdot)); 298 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 299 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 300 } 301 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 302 PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 303 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 304 PetscFunctionReturn(0); 305 } 306 307 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems *PetscOptionsObject) 308 { 309 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 310 PetscBool flg = PETSC_FALSE; 311 PetscViewer viewer; 312 313 PetscFunctionBegin; 314 PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 315 PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 316 if (flg) { 317 PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 318 PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 319 } 320 flg = pseudo->increment_dt_from_initial_dt; 321 PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 322 pseudo->increment_dt_from_initial_dt = flg; 323 PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 324 PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 325 PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 326 PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 327 PetscOptionsHeadEnd(); 328 PetscFunctionReturn(0); 329 } 330 331 static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 332 { 333 PetscBool isascii; 334 335 PetscFunctionBegin; 336 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 337 if (isascii) { 338 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 339 PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 340 PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 341 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 342 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 343 PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 344 } 345 PetscFunctionReturn(0); 346 } 347 348 /* ----------------------------------------------------------------------------- */ 349 /*@C 350 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 351 last timestep. 352 353 Logically Collective on ts 354 355 Input Parameters: 356 + ts - timestep context 357 . dt - user-defined function to verify timestep 358 - ctx - [optional] user-defined context for private data 359 for the timestep verification routine (may be NULL) 360 361 Calling sequence of func: 362 $ func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 363 364 + update - latest solution vector 365 . ctx - [optional] timestep context 366 . newdt - the timestep to use for the next step 367 - flag - flag indicating whether the last time step was acceptable 368 369 Level: advanced 370 371 Note: 372 The routine set here will be called by `TSPseudoVerifyTimeStep()` 373 during the timestepping process. 374 375 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 376 @*/ 377 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS, Vec, void *, PetscReal *, PetscBool *), void *ctx) 378 { 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 381 PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode(*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 382 PetscFunctionReturn(0); 383 } 384 385 /*@ 386 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 387 dt when using the TSPseudoTimeStepDefault() routine. 388 389 Logically Collective on ts 390 391 Input Parameters: 392 + ts - the timestep context 393 - inc - the scaling factor >= 1.0 394 395 Options Database Key: 396 . -ts_pseudo_increment <increment> - set pseudo increment 397 398 Level: advanced 399 400 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 401 @*/ 402 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 403 { 404 PetscFunctionBegin; 405 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 406 PetscValidLogicalCollectiveReal(ts, inc, 2); 407 PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 408 PetscFunctionReturn(0); 409 } 410 411 /*@ 412 TSPseudoSetMaxTimeStep - Sets the maximum time step 413 when using the TSPseudoTimeStepDefault() routine. 414 415 Logically Collective on ts 416 417 Input Parameters: 418 + ts - the timestep context 419 - maxdt - the maximum time step, use a non-positive value to deactivate 420 421 Options Database Key: 422 . -ts_pseudo_max_dt <increment> - set pseudo max dt 423 424 Level: advanced 425 426 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 427 @*/ 428 PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 429 { 430 PetscFunctionBegin; 431 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 432 PetscValidLogicalCollectiveReal(ts, maxdt, 2); 433 PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 434 PetscFunctionReturn(0); 435 } 436 437 /*@ 438 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 439 is computed via the formula 440 $ dt = initial_dt*initial_fnorm/current_fnorm 441 rather than the default update, 442 $ dt = current_dt*previous_fnorm/current_fnorm. 443 444 Logically Collective on ts 445 446 Input Parameter: 447 . ts - the timestep context 448 449 Options Database Key: 450 . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment 451 452 Level: advanced 453 454 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 455 @*/ 456 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 457 { 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 460 PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 461 PetscFunctionReturn(0); 462 } 463 464 /*@C 465 TSPseudoSetTimeStep - Sets the user-defined routine to be 466 called at each pseudo-timestep to update the timestep. 467 468 Logically Collective on ts 469 470 Input Parameters: 471 + ts - timestep context 472 . dt - function to compute timestep 473 - ctx - [optional] user-defined context for private data 474 required by the function (may be NULL) 475 476 Calling sequence of func: 477 $ func (TS ts,PetscReal *newdt,void *ctx); 478 479 + newdt - the newly computed timestep 480 - ctx - [optional] timestep context 481 482 Level: intermediate 483 484 Notes: 485 The routine set here will be called by `TSPseudoComputeTimeStep()` 486 during the timestepping process. 487 488 If not set then `TSPseudoTimeStepDefault()` is automatically used 489 490 .seealso: [](chapter_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 491 @*/ 492 PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS, PetscReal *, void *), void *ctx) 493 { 494 PetscFunctionBegin; 495 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 496 PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode(*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 497 PetscFunctionReturn(0); 498 } 499 500 /* ----------------------------------------------------------------------------- */ 501 502 typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 503 static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 504 { 505 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 506 507 PetscFunctionBegin; 508 pseudo->verify = dt; 509 pseudo->verifyctx = ctx; 510 PetscFunctionReturn(0); 511 } 512 513 static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 514 { 515 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 516 517 PetscFunctionBegin; 518 pseudo->dt_increment = inc; 519 PetscFunctionReturn(0); 520 } 521 522 static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 523 { 524 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 525 526 PetscFunctionBegin; 527 pseudo->dt_max = maxdt; 528 PetscFunctionReturn(0); 529 } 530 531 static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 532 { 533 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 534 535 PetscFunctionBegin; 536 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 537 PetscFunctionReturn(0); 538 } 539 540 typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 541 static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 542 { 543 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 544 545 PetscFunctionBegin; 546 pseudo->dt = dt; 547 pseudo->dtctx = ctx; 548 PetscFunctionReturn(0); 549 } 550 551 /* ----------------------------------------------------------------------------- */ 552 /*MC 553 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 554 555 This method solves equations of the form 556 557 $ F(X,Xdot) = 0 558 559 for steady state using the iteration 560 561 $ [G'] S = -F(X,0) 562 $ X += S 563 564 where 565 566 $ G(Y) = F(Y,(Y-X)/dt) 567 568 This is linearly-implicit Euler with the residual always evaluated "at steady 569 state". See note below. 570 571 Options Database Keys: 572 + -ts_pseudo_increment <real> - ratio of increase dt 573 . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 574 . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 575 - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 576 577 Level: beginner 578 579 Notes: 580 The residual computed by this method includes the transient term (Xdot is computed instead of 581 always being zero), but since the prediction from the last step is always the solution from the 582 last step, on the first Newton iteration we have 583 584 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 585 586 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 587 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 588 algorithm is no longer the one described in the referenced papers. 589 590 References: 591 + * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 592 - * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 593 594 .seealso: [](chapter_ts), `TSCreate()`, `TS`, `TSSetType()` 595 M*/ 596 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 597 { 598 TS_Pseudo *pseudo; 599 SNES snes; 600 SNESType stype; 601 602 PetscFunctionBegin; 603 ts->ops->reset = TSReset_Pseudo; 604 ts->ops->destroy = TSDestroy_Pseudo; 605 ts->ops->view = TSView_Pseudo; 606 ts->ops->setup = TSSetUp_Pseudo; 607 ts->ops->step = TSStep_Pseudo; 608 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 609 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 610 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 611 ts->default_adapt_type = TSADAPTNONE; 612 ts->usessnes = PETSC_TRUE; 613 614 PetscCall(TSGetSNES(ts, &snes)); 615 PetscCall(SNESGetType(snes, &stype)); 616 if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 617 618 PetscCall(PetscNew(&pseudo)); 619 ts->data = (void *)pseudo; 620 621 pseudo->dt = TSPseudoTimeStepDefault; 622 pseudo->dtctx = NULL; 623 pseudo->dt_increment = 1.1; 624 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 625 pseudo->fnorm = -1; 626 pseudo->fnorm_initial = -1; 627 pseudo->fnorm_previous = -1; 628 #if defined(PETSC_USE_REAL_SINGLE) 629 pseudo->fatol = 1.e-25; 630 pseudo->frtol = 1.e-5; 631 #else 632 pseudo->fatol = 1.e-50; 633 pseudo->frtol = 1.e-12; 634 #endif 635 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 636 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 637 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 638 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 639 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 640 PetscFunctionReturn(0); 641 } 642 643 /*@C 644 TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 645 646 Collective on ts 647 648 Input Parameters: 649 + ts - the timestep context 650 - dtctx - unused timestep context 651 652 Output Parameter: 653 . newdt - the timestep to use for the next step 654 655 Level: advanced 656 657 .seealso: [](chapter_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 658 @*/ 659 PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 660 { 661 TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 662 PetscReal inc = pseudo->dt_increment; 663 664 PetscFunctionBegin; 665 PetscCall(VecZeroEntries(pseudo->xdot)); 666 PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 667 PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 668 if (pseudo->fnorm_initial < 0) { 669 /* first time through so compute initial function norm */ 670 pseudo->fnorm_initial = pseudo->fnorm; 671 pseudo->fnorm_previous = pseudo->fnorm; 672 } 673 if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 674 else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 675 else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 676 if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 677 pseudo->fnorm_previous = pseudo->fnorm; 678 PetscFunctionReturn(0); 679 } 680