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