xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 940ebdd4d53c9bd2466ad09f7a10e4b0eba2bd47)
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 /*@
408   TSPseudoComputeTimeStep - Computes the next timestep for a currently running
409   pseudo-timestepping process.
410 
411   Collective
412 
413   Input Parameter:
414 . ts - timestep context
415 
416   Output Parameter:
417 . dt - newly computed timestep
418 
419   Level: developer
420 
421   Note:
422   The routine to be called here to compute the timestep should be
423   set by calling `TSPseudoSetTimeStep()`.
424 
425 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
426 @*/
427 PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt)
428 {
429   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
430 
431   PetscFunctionBegin;
432   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
433   PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx));
434   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
435   PetscFunctionReturn(PETSC_SUCCESS);
436 }
437 
438 /*@C
439   TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
440 
441   Collective, No Fortran Support
442 
443   Input Parameters:
444 + ts     - the timestep context
445 . dtctx  - unused timestep context
446 - update - latest solution vector
447 
448   Output Parameters:
449 + newdt - the timestep to use for the next step
450 - flag  - flag indicating whether the last time step was acceptable
451 
452   Level: advanced
453 
454   Note:
455   This routine always returns a flag of 1, indicating an acceptable
456   timestep.
457 
458 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
459 @*/
460 PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag)
461 {
462   PetscFunctionBegin;
463   // NOTE: This function was never used
464   *flag = PETSC_TRUE;
465   PetscFunctionReturn(PETSC_SUCCESS);
466 }
467 
468 /*@
469   TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
470 
471   Collective
472 
473   Input Parameters:
474 + ts     - timestep context
475 - update - latest solution vector
476 
477   Output Parameters:
478 + dt   - newly computed timestep (if it had to shrink)
479 - flag - indicates if current timestep was ok
480 
481   Level: advanced
482 
483   Notes:
484   The routine to be called here to compute the timestep should be
485   set by calling `TSPseudoSetVerifyTimeStep()`.
486 
487 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
488 @*/
489 PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag)
490 {
491   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
492 
493   PetscFunctionBegin;
494   // NOTE: This function is never used
495   *flag = PETSC_TRUE;
496   if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag));
497   PetscFunctionReturn(PETSC_SUCCESS);
498 }
499 
500 /*@C
501   TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
502   last timestep.
503 
504   Logically Collective
505 
506   Input Parameters:
507 + ts  - timestep context
508 . dt  - user-defined function to verify timestep
509 - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`)
510 
511   Calling sequence of `func`:
512 + ts     - the time-step context
513 . update - latest solution vector
514 . ctx    - [optional] user-defined timestep context
515 . newdt  - the timestep to use for the next step
516 - flag   - flag indicating whether the last time step was acceptable
517 
518   Level: advanced
519 
520   Note:
521   The routine set here will be called by `TSPseudoVerifyTimeStep()`
522   during the timestepping process.
523 
524 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
525 @*/
526 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx)
527 {
528   PetscFunctionBegin;
529   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
530   PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx));
531   PetscFunctionReturn(PETSC_SUCCESS);
532 }
533 
534 /*@
535   TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
536   dt when using the TSPseudoTimeStepDefault() routine.
537 
538   Logically Collective
539 
540   Input Parameters:
541 + ts  - the timestep context
542 - inc - the scaling factor >= 1.0
543 
544   Options Database Key:
545 . -ts_pseudo_increment <increment> - set pseudo increment
546 
547   Level: advanced
548 
549 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
550 @*/
551 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc)
552 {
553   PetscFunctionBegin;
554   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
555   PetscValidLogicalCollectiveReal(ts, inc, 2);
556   PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc));
557   PetscFunctionReturn(PETSC_SUCCESS);
558 }
559 
560 /*@
561   TSPseudoSetMaxTimeStep - Sets the maximum time step
562   when using the TSPseudoTimeStepDefault() routine.
563 
564   Logically Collective
565 
566   Input Parameters:
567 + ts    - the timestep context
568 - maxdt - the maximum time step, use a non-positive value to deactivate
569 
570   Options Database Key:
571 . -ts_pseudo_max_dt <increment> - set pseudo max dt
572 
573   Level: advanced
574 
575 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
576 @*/
577 PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt)
578 {
579   PetscFunctionBegin;
580   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
581   PetscValidLogicalCollectiveReal(ts, maxdt, 2);
582   PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
583   PetscFunctionReturn(PETSC_SUCCESS);
584 }
585 
586 /*@
587   TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
588   is computed via the formula $  dt = initial\_dt*initial\_fnorm/current\_fnorm $  rather than the default update, $  dt = current\_dt*previous\_fnorm/current\_fnorm.$
589 
590   Logically Collective
591 
592   Input Parameter:
593 . ts - the timestep context
594 
595   Options Database Key:
596 . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment
597 
598   Level: advanced
599 
600 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
601 @*/
602 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
603 {
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
606   PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts));
607   PetscFunctionReturn(PETSC_SUCCESS);
608 }
609 
610 /*@C
611   TSPseudoSetTimeStep - Sets the user-defined routine to be
612   called at each pseudo-timestep to update the timestep.
613 
614   Logically Collective
615 
616   Input Parameters:
617 + ts  - timestep context
618 . dt  - function to compute timestep
619 - ctx - [optional] user-defined context for private data required by the function (may be `NULL`)
620 
621   Calling sequence of `dt`:
622 + ts    - the `TS` context
623 . newdt - the newly computed timestep
624 - ctx   - [optional] user-defined context
625 
626   Level: intermediate
627 
628   Notes:
629   The routine set here will be called by `TSPseudoComputeTimeStep()`
630   during the timestepping process.
631 
632   If not set then `TSPseudoTimeStepDefault()` is automatically used
633 
634 .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
635 @*/
636 PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx)
637 {
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
640   PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx));
641   PetscFunctionReturn(PETSC_SUCCESS);
642 }
643 
644 typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/
645 static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx)
646 {
647   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
648 
649   PetscFunctionBegin;
650   pseudo->verify    = dt;
651   pseudo->verifyctx = ctx;
652   PetscFunctionReturn(PETSC_SUCCESS);
653 }
654 
655 static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc)
656 {
657   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
658 
659   PetscFunctionBegin;
660   pseudo->dt_increment = inc;
661   PetscFunctionReturn(PETSC_SUCCESS);
662 }
663 
664 static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt)
665 {
666   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
667 
668   PetscFunctionBegin;
669   pseudo->dt_max = maxdt;
670   PetscFunctionReturn(PETSC_SUCCESS);
671 }
672 
673 static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
674 {
675   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
676 
677   PetscFunctionBegin;
678   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
679   PetscFunctionReturn(PETSC_SUCCESS);
680 }
681 
682 typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/
683 static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx)
684 {
685   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
686 
687   PetscFunctionBegin;
688   pseudo->dt    = dt;
689   pseudo->dtctx = ctx;
690   PetscFunctionReturn(PETSC_SUCCESS);
691 }
692 
693 /*MC
694   TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97`
695 
696   This method solves equations of the form
697 
698   $$
699   F(X,Xdot) = 0
700   $$
701 
702   for steady state using the iteration
703 
704   $$
705   [G'] S = -F(X,0)
706   X += S
707   $$
708 
709   where
710 
711   $$
712   G(Y) = F(Y,(Y-X)/dt)
713   $$
714 
715   This is linearly-implicit Euler with the residual always evaluated "at steady
716   state".  See note below.
717 
718   In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`.
719   It determines the next timestep via
720 
721   $$
722   dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert}
723   $$
724 
725   where $r$ is an additional growth factor (set by `-ts_pseudo_increment`).
726   An alternative formulation is also available that uses the initial timestep and function norm.
727 
728   $$
729   dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert}
730   $$
731 
732   This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`.
733   For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`.
734 
735   Options Database Keys:
736 +  -ts_pseudo_increment <real>                     - ratio of increase dt
737 .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
738 .  -ts_pseudo_max_dt                               - Maximum dt for adaptive timestepping algorithm
739 .  -ts_pseudo_monitor                              - Monitor convergence of the function norm
740 .  -ts_pseudo_fatol <atol>                         - stop iterating when the function norm is less than `atol`
741 -  -ts_pseudo_frtol <rtol>                         - stop iterating when the function norm divided by the initial function norm is less than `rtol`
742 
743   Level: beginner
744 
745   Notes:
746   The residual computed by this method includes the transient term (Xdot is computed instead of
747   always being zero), but since the prediction from the last step is always the solution from the
748   last step, on the first Newton iteration we have
749 
750   $$
751   Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
752   $$
753 
754   The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it
755   is $ dF/dX + shift*1/dt $. Hence  still contains the $ 1/dt $ term so the Jacobian is not simply the
756   Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$.
757 
758   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
759   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
760   algorithm is no longer the one described in the referenced papers.
761   By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers.
762   Setting the `SNESType` via `-snes_type` will override this default setting.
763 
764 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`
765 M*/
766 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
767 {
768   TS_Pseudo *pseudo;
769   SNES       snes;
770   SNESType   stype;
771 
772   PetscFunctionBegin;
773   ts->ops->reset          = TSReset_Pseudo;
774   ts->ops->destroy        = TSDestroy_Pseudo;
775   ts->ops->view           = TSView_Pseudo;
776   ts->ops->setup          = TSSetUp_Pseudo;
777   ts->ops->step           = TSStep_Pseudo;
778   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
779   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
780   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
781   ts->default_adapt_type  = TSADAPTTSPSEUDO;
782   ts->usessnes            = PETSC_TRUE;
783 
784   PetscCall(TSAdaptRegister(TSADAPTTSPSEUDO, TSAdaptCreate_TSPseudo));
785 
786   PetscCall(TSGetSNES(ts, &snes));
787   PetscCall(SNESGetType(snes, &stype));
788   if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY));
789 
790   PetscCall(PetscNew(&pseudo));
791   ts->data = (void *)pseudo;
792 
793   pseudo->dt                           = TSPseudoTimeStepDefault;
794   pseudo->dtctx                        = NULL;
795   pseudo->dt_increment                 = 1.1;
796   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
797   pseudo->fnorm                        = -1;
798   pseudo->fnorm_initial                = -1;
799   pseudo->fnorm_previous               = -1;
800 #if defined(PETSC_USE_REAL_SINGLE)
801   pseudo->fatol = 1.e-25;
802   pseudo->frtol = 1.e-5;
803 #else
804   pseudo->fatol = 1.e-50;
805   pseudo->frtol = 1.e-12;
806 #endif
807   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo));
808   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo));
809   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo));
810   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo));
811   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo));
812   PetscFunctionReturn(PETSC_SUCCESS);
813 }
814