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