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