xref: /petsc/src/ts/impls/pseudo/posindep.c (revision eb636060ec67a08151583c50986cf2393a5eae68)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
52d3f70b5SBarry Smith 
62d3f70b5SBarry Smith typedef struct {
72d3f70b5SBarry Smith   Vec update; /* work vector where new solution is formed */
82d3f70b5SBarry Smith   Vec func;   /* work vector where F(t[i],u[i]) is stored */
96f2d6a7bSJed Brown   Vec xdot;   /* work vector for time derivative of state */
102d3f70b5SBarry Smith 
112d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
122d3f70b5SBarry Smith 
136849ba73SBarry Smith   PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */
142d3f70b5SBarry Smith   void *dtctx;
15ace3abfcSBarry Smith   PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */
167bf11e45SBarry Smith   void *verifyctx;
172d3f70b5SBarry Smith 
18cdbf8f93SLisandro Dalcin   PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */
1987828ca2SBarry Smith   PetscReal fnorm_previous;
2028aa8177SBarry Smith 
21cdbf8f93SLisandro Dalcin   PetscReal dt_initial;   /* initial time-step */
2287828ca2SBarry Smith   PetscReal dt_increment; /* scaling that dt is incremented each time-step */
2386534af1SJed Brown   PetscReal dt_max;       /* maximum time step */
24ace3abfcSBarry Smith   PetscBool increment_dt_from_initial_dt;
253118ae5eSBarry Smith   PetscReal fatol, frtol;
26*eb636060SBarry Smith 
27*eb636060SBarry Smith   PetscObjectState Xstate; /* state of input vector to TSComputeIFunction() with 0 Xdot */
287bf11e45SBarry Smith } TS_Pseudo;
292d3f70b5SBarry Smith 
302d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
312d3f70b5SBarry Smith 
32cc4c1da9SBarry Smith /*@
337bf11e45SBarry Smith   TSPseudoComputeTimeStep - Computes the next timestep for a currently running
34564e8f4eSLois Curfman McInnes   pseudo-timestepping process.
352d3f70b5SBarry Smith 
36c3339decSBarry Smith   Collective
3715091d37SBarry Smith 
387bf11e45SBarry Smith   Input Parameter:
397bf11e45SBarry Smith . ts - timestep context
407bf11e45SBarry Smith 
417bf11e45SBarry Smith   Output Parameter:
42fb4a63b6SLois Curfman McInnes . dt - newly computed timestep
43fb4a63b6SLois Curfman McInnes 
448d359177SBarry Smith   Level: developer
45564e8f4eSLois Curfman McInnes 
46bcf0153eSBarry Smith   Note:
47564e8f4eSLois Curfman McInnes   The routine to be called here to compute the timestep should be
48bcf0153eSBarry Smith   set by calling `TSPseudoSetTimeStep()`.
49564e8f4eSLois Curfman McInnes 
501cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
517bf11e45SBarry Smith @*/
52d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt)
53d71ae5a4SJacob Faibussowitsch {
547bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
557bf11e45SBarry Smith 
563a40ed3dSBarry Smith   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
589566063dSJacob Faibussowitsch   PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx));
599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617bf11e45SBarry Smith }
627bf11e45SBarry Smith 
637bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
647bf11e45SBarry Smith /*@C
658d359177SBarry Smith   TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
667bf11e45SBarry Smith 
67cc4c1da9SBarry Smith   Collective, No Fortran Support
6815091d37SBarry Smith 
697bf11e45SBarry Smith   Input Parameters:
7015091d37SBarry Smith + ts     - the timestep context
717bf11e45SBarry Smith . dtctx  - unused timestep context
7215091d37SBarry Smith - update - latest solution vector
737bf11e45SBarry Smith 
74564e8f4eSLois Curfman McInnes   Output Parameters:
7515091d37SBarry Smith + newdt - the timestep to use for the next step
7615091d37SBarry Smith - flag  - flag indicating whether the last time step was acceptable
777bf11e45SBarry Smith 
7815091d37SBarry Smith   Level: advanced
79fee21e36SBarry Smith 
80564e8f4eSLois Curfman McInnes   Note:
81564e8f4eSLois Curfman McInnes   This routine always returns a flag of 1, indicating an acceptable
82564e8f4eSLois Curfman McInnes   timestep.
83564e8f4eSLois Curfman McInnes 
841cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
857bf11e45SBarry Smith @*/
86d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag)
87d71ae5a4SJacob Faibussowitsch {
883a40ed3dSBarry Smith   PetscFunctionBegin;
89a7cc72afSBarry Smith   *flag = PETSC_TRUE;
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
917bf11e45SBarry Smith }
927bf11e45SBarry Smith 
937bf11e45SBarry Smith /*@
94564e8f4eSLois Curfman McInnes   TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
957bf11e45SBarry Smith 
96c3339decSBarry Smith   Collective
9715091d37SBarry Smith 
98fb4a63b6SLois Curfman McInnes   Input Parameters:
9915091d37SBarry Smith + ts     - timestep context
10015091d37SBarry Smith - update - latest solution vector
1017bf11e45SBarry Smith 
102fb4a63b6SLois Curfman McInnes   Output Parameters:
10315091d37SBarry Smith + dt   - newly computed timestep (if it had to shrink)
10415091d37SBarry Smith - flag - indicates if current timestep was ok
1057bf11e45SBarry Smith 
10615091d37SBarry Smith   Level: advanced
107fee21e36SBarry Smith 
108564e8f4eSLois Curfman McInnes   Notes:
109564e8f4eSLois Curfman McInnes   The routine to be called here to compute the timestep should be
110bcf0153eSBarry Smith   set by calling `TSPseudoSetVerifyTimeStep()`.
111564e8f4eSLois Curfman McInnes 
1121cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
1137bf11e45SBarry Smith @*/
114d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag)
115d71ae5a4SJacob Faibussowitsch {
1167bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1177bf11e45SBarry Smith 
1183a40ed3dSBarry Smith   PetscFunctionBegin;
119cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
1201baa6e33SBarry Smith   if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1227bf11e45SBarry Smith }
1237bf11e45SBarry Smith 
1247bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1257bf11e45SBarry Smith 
126d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts)
127d71ae5a4SJacob Faibussowitsch {
128277b19d0SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
129be5899b3SLisandro Dalcin   PetscInt   nits, lits, reject;
130cdbf8f93SLisandro Dalcin   PetscBool  stepok;
131be5899b3SLisandro Dalcin   PetscReal  next_time_step = ts->time_step;
1322d3f70b5SBarry Smith 
1333a40ed3dSBarry Smith   PetscFunctionBegin;
134*eb636060SBarry Smith   if (ts->steps == 0) {
135*eb636060SBarry Smith     pseudo->dt_initial = ts->time_step;
136*eb636060SBarry Smith     PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* in all future updates pseudo->update already contains the current time solution */
137*eb636060SBarry Smith   }
1389566063dSJacob Faibussowitsch   PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step));
139cdbf8f93SLisandro Dalcin   for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) {
140cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
1419566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime + ts->time_step));
1429566063dSJacob Faibussowitsch     PetscCall(SNESSolve(ts->snes, NULL, pseudo->update));
1439566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1449566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1459371c9d4SSatish Balay     ts->snes_its += nits;
1469371c9d4SSatish Balay     ts->ksp_its += lits;
147f4f49eeaSPierre Jolivet     PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &pseudo->update));
1489566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok));
1499371c9d4SSatish Balay     if (!stepok) {
1509371c9d4SSatish Balay       next_time_step = ts->time_step;
1519371c9d4SSatish Balay       continue;
1529371c9d4SSatish Balay     }
153*eb636060SBarry Smith     pseudo->fnorm = -1; /* The current norm is no longer valid */
1549566063dSJacob Faibussowitsch     PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok));
155cdbf8f93SLisandro Dalcin     if (stepok) break;
156cdbf8f93SLisandro Dalcin   }
157be5899b3SLisandro Dalcin   if (reject >= ts->max_reject) {
158be5899b3SLisandro Dalcin     ts->reason = TS_DIVERGED_STEP_REJECTED;
15963a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject));
1603ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1617bf11e45SBarry Smith   }
162be5899b3SLisandro Dalcin 
1639566063dSJacob Faibussowitsch   PetscCall(VecCopy(pseudo->update, ts->vec_sol));
164be5899b3SLisandro Dalcin   ts->ptime += ts->time_step;
165be5899b3SLisandro Dalcin   ts->time_step = next_time_step;
166be5899b3SLisandro Dalcin 
1679566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(pseudo->xdot));
168*eb636060SBarry Smith   PetscCall(TSComputeIFunction(ts, ts->ptime, pseudo->update, pseudo->xdot, pseudo->func, PETSC_FALSE));
169*eb636060SBarry Smith   PetscCall(PetscObjectStateGet((PetscObject)pseudo->update, &pseudo->Xstate));
1709566063dSJacob Faibussowitsch   PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
171*eb636060SBarry Smith 
1723118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1733118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
17463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol));
1753ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1763118ae5eSBarry Smith   }
1773118ae5eSBarry Smith   if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) {
1783118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
17963a3b9bcSJacob Faibussowitsch     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));
1803ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1813118ae5eSBarry Smith   }
1823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1832d3f70b5SBarry Smith }
1842d3f70b5SBarry Smith 
1852d3f70b5SBarry Smith /*------------------------------------------------------------*/
186d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts)
187d71ae5a4SJacob Faibussowitsch {
1887bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1892d3f70b5SBarry Smith 
1903a40ed3dSBarry Smith   PetscFunctionBegin;
1919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->update));
1929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->func));
1939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->xdot));
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1952d3f70b5SBarry Smith }
1962d3f70b5SBarry Smith 
197d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts)
198d71ae5a4SJacob Faibussowitsch {
199277b19d0SLisandro Dalcin   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(TSReset_Pseudo(ts));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
2029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL));
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL));
2059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL));
2069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208277b19d0SLisandro Dalcin }
2092d3f70b5SBarry Smith 
2102d3f70b5SBarry Smith /*------------------------------------------------------------*/
2112d3f70b5SBarry Smith 
2126f2d6a7bSJed Brown /*
2136f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2146f2d6a7bSJed Brown */
215d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot)
216d71ae5a4SJacob Faibussowitsch {
2176f2d6a7bSJed Brown   TS_Pseudo        *pseudo = (TS_Pseudo *)ts->data;
218193ac0bcSJed Brown   const PetscScalar mdt    = 1.0 / ts->time_step, *xnp1, *xn;
219193ac0bcSJed Brown   PetscScalar      *xdot;
220a7cc72afSBarry Smith   PetscInt          i, n;
2212d3f70b5SBarry Smith 
2223a40ed3dSBarry Smith   PetscFunctionBegin;
223aab5bcd8SJed Brown   *Xdot = NULL;
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ts->vec_sol, &xn));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xnp1));
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pseudo->xdot, &xdot));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
228bbd56ea5SKarl Rupp   for (i = 0; i < n; i++) xdot[i] = mdt * (xnp1[i] - xn[i]);
2299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ts->vec_sol, &xn));
2309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xnp1));
2319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pseudo->xdot, &xdot));
2326f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2342d3f70b5SBarry Smith }
2352d3f70b5SBarry Smith 
2366f2d6a7bSJed Brown /*
2376f2d6a7bSJed Brown     The transient residual is
2386f2d6a7bSJed Brown 
2396f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2406f2d6a7bSJed Brown 
2416f2d6a7bSJed Brown     or for ODE,
2426f2d6a7bSJed Brown 
2436f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2446f2d6a7bSJed Brown 
2456f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2466f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2476f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2486f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2496f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
250*eb636060SBarry Smith     algorithm, it only takes this one full Newton step with the steady state
2516f2d6a7bSJed Brown     residual, and then advances to the next time step.
252*eb636060SBarry Smith 
253*eb636060SBarry Smith     This routine avoids a repeated identical call to TSComputeRHSFunction() when that result
254*eb636060SBarry Smith     is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo()
255*eb636060SBarry Smith 
2566f2d6a7bSJed Brown */
257d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts)
258d71ae5a4SJacob Faibussowitsch {
2596f2d6a7bSJed Brown   Vec              Xdot;
260*eb636060SBarry Smith   TSIFunctionFn   *ifunction;
261*eb636060SBarry Smith   TSRHSFunctionFn *rhsfunction;
262*eb636060SBarry Smith   void            *ctx;
263*eb636060SBarry Smith   DM               dm;
264*eb636060SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo *)ts->data;
265*eb636060SBarry Smith   PetscBool        KSPSNES;
266*eb636060SBarry Smith   PetscObjectState Xstate;
2672d3f70b5SBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionBegin;
269*eb636060SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
270*eb636060SBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
271*eb636060SBarry Smith   PetscValidHeaderSpecific(Y, VEC_CLASSID, 3);
272*eb636060SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 4);
273*eb636060SBarry Smith 
274*eb636060SBarry Smith   PetscCall(TSGetDM(ts, &dm));
275*eb636060SBarry Smith   PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx));
276*eb636060SBarry Smith   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
277*eb636060SBarry Smith   PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()");
278*eb636060SBarry Smith 
279*eb636060SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESKSPONLY, &KSPSNES));
2809566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
281*eb636060SBarry Smith   PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate));
282*eb636060SBarry Smith   if (Xstate != pseudo->Xstate || ifunction || !KSPSNES) {
2839566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, Xdot, Y, PETSC_FALSE));
284*eb636060SBarry Smith   } else {
285*eb636060SBarry Smith     /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */
286*eb636060SBarry Smith     /* note that pseudo->func contains the negation of TSComputeRHSFunction() */
287*eb636060SBarry Smith     PetscCall(VecWAXPY(Y, 1, pseudo->func, Xdot));
288*eb636060SBarry Smith   }
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2906f2d6a7bSJed Brown }
2912d3f70b5SBarry Smith 
2926f2d6a7bSJed Brown /*
2936f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2946f2d6a7bSJed Brown 
2956f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2966f2d6a7bSJed Brown 
2976f2d6a7bSJed Brown     and for ODE:
2986f2d6a7bSJed Brown 
2996f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
3006f2d6a7bSJed Brown */
301d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts)
302d71ae5a4SJacob Faibussowitsch {
3036f2d6a7bSJed Brown   Vec Xdot;
3046f2d6a7bSJed Brown 
3056f2d6a7bSJed Brown   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
3079566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE));
3083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3092d3f70b5SBarry Smith }
3102d3f70b5SBarry Smith 
311d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts)
312d71ae5a4SJacob Faibussowitsch {
3137bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
3142d3f70b5SBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update));
3179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func));
3189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot));
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3202d3f70b5SBarry Smith }
3212d3f70b5SBarry Smith /*------------------------------------------------------------*/
3222d3f70b5SBarry Smith 
323d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy)
324d71ae5a4SJacob Faibussowitsch {
3257bf11e45SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
326ce94432eSBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
3272d3f70b5SBarry Smith 
3283a40ed3dSBarry Smith   PetscFunctionBegin;
329193ac0bcSJed Brown   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
3309566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
3319566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
332*eb636060SBarry Smith     PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate));
3339566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
334193ac0bcSJed Brown   }
3359566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
33663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm));
3379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3392d3f70b5SBarry Smith }
3402d3f70b5SBarry Smith 
341ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject)
342d71ae5a4SJacob Faibussowitsch {
3434bbc92c1SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
344ace3abfcSBarry Smith   PetscBool   flg    = PETSC_FALSE;
345649052a6SBarry Smith   PetscViewer viewer;
3462d3f70b5SBarry Smith 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
348d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options");
3499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL));
3502d3f70b5SBarry Smith   if (flg) {
3519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer));
35249abdd8aSBarry Smith     PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
35328aa8177SBarry Smith   }
354be5899b3SLisandro Dalcin   flg = pseudo->increment_dt_from_initial_dt;
3559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL));
356be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
3579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL));
3589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL));
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL));
361d0609cedSBarry Smith   PetscOptionsHeadEnd();
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3632d3f70b5SBarry Smith }
3642d3f70b5SBarry Smith 
365d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer)
366d71ae5a4SJacob Faibussowitsch {
3673118ae5eSBarry Smith   PetscBool isascii;
368d52bd9f3SBarry Smith 
3693a40ed3dSBarry Smith   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3713118ae5eSBarry Smith   if (isascii) {
3723118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
3739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  frtol - relative tolerance in function value %g\n", (double)pseudo->frtol));
3749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol));
3759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_initial - initial timestep %g\n", (double)pseudo->dt_initial));
3769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment));
3779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_max - maximum time %g\n", (double)pseudo->dt_max));
3783118ae5eSBarry Smith   }
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3802d3f70b5SBarry Smith }
3812d3f70b5SBarry Smith 
382ac226902SBarry Smith /*@C
38382bf6240SBarry Smith   TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
38482bf6240SBarry Smith   last timestep.
38582bf6240SBarry Smith 
386c3339decSBarry Smith   Logically Collective
38715091d37SBarry Smith 
38882bf6240SBarry Smith   Input Parameters:
38915091d37SBarry Smith + ts  - timestep context
39082bf6240SBarry Smith . dt  - user-defined function to verify timestep
3912fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`)
39282bf6240SBarry Smith 
39320f4b53cSBarry Smith   Calling sequence of `func`:
3942fe279fdSBarry Smith + ts     - the time-step context
3952fe279fdSBarry Smith . update - latest solution vector
39614d0ab18SJacob Faibussowitsch . ctx    - [optional] user-defined timestep context
39782bf6240SBarry Smith . newdt  - the timestep to use for the next step
398a2b725a8SWilliam Gropp - flag   - flag indicating whether the last time step was acceptable
39982bf6240SBarry Smith 
400bcf0153eSBarry Smith   Level: advanced
401bcf0153eSBarry Smith 
402bcf0153eSBarry Smith   Note:
403bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoVerifyTimeStep()`
40482bf6240SBarry Smith   during the timestepping process.
40582bf6240SBarry Smith 
4061cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
40782bf6240SBarry Smith @*/
40814d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx)
409d71ae5a4SJacob Faibussowitsch {
41082bf6240SBarry Smith   PetscFunctionBegin;
4110700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
412cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx));
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41482bf6240SBarry Smith }
41582bf6240SBarry Smith 
41682bf6240SBarry Smith /*@
41782bf6240SBarry Smith   TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
4188d359177SBarry Smith   dt when using the TSPseudoTimeStepDefault() routine.
41982bf6240SBarry Smith 
420c3339decSBarry Smith   Logically Collective
421fee21e36SBarry Smith 
42215091d37SBarry Smith   Input Parameters:
42315091d37SBarry Smith + ts  - the timestep context
42415091d37SBarry Smith - inc - the scaling factor >= 1.0
42515091d37SBarry Smith 
42682bf6240SBarry Smith   Options Database Key:
42767b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment
42882bf6240SBarry Smith 
42915091d37SBarry Smith   Level: advanced
43015091d37SBarry Smith 
4311cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
43282bf6240SBarry Smith @*/
433d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc)
434d71ae5a4SJacob Faibussowitsch {
43582bf6240SBarry Smith   PetscFunctionBegin;
4360700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
437c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts, inc, 2);
438cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc));
4393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44082bf6240SBarry Smith }
44182bf6240SBarry Smith 
44286534af1SJed Brown /*@
44386534af1SJed Brown   TSPseudoSetMaxTimeStep - Sets the maximum time step
4448d359177SBarry Smith   when using the TSPseudoTimeStepDefault() routine.
44586534af1SJed Brown 
446c3339decSBarry Smith   Logically Collective
44786534af1SJed Brown 
44886534af1SJed Brown   Input Parameters:
44986534af1SJed Brown + ts    - the timestep context
45086534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate
45186534af1SJed Brown 
45286534af1SJed Brown   Options Database Key:
45367b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt
45486534af1SJed Brown 
45586534af1SJed Brown   Level: advanced
45686534af1SJed Brown 
4571cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
45886534af1SJed Brown @*/
459d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt)
460d71ae5a4SJacob Faibussowitsch {
46186534af1SJed Brown   PetscFunctionBegin;
46286534af1SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
46386534af1SJed Brown   PetscValidLogicalCollectiveReal(ts, maxdt, 2);
464cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46686534af1SJed Brown }
46786534af1SJed Brown 
46882bf6240SBarry Smith /*@
46982bf6240SBarry Smith   TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
470b44f4de4SBarry Smith   is computed via the formula $  dt = initial\_dt*initial\_fnorm/current\_fnorm $  rather than the default update, $  dt = current\_dt*previous\_fnorm/current\_fnorm.$
47182bf6240SBarry Smith 
472c3339decSBarry Smith   Logically Collective
47315091d37SBarry Smith 
47482bf6240SBarry Smith   Input Parameter:
47582bf6240SBarry Smith . ts - the timestep context
47682bf6240SBarry Smith 
47782bf6240SBarry Smith   Options Database Key:
478b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment
47982bf6240SBarry Smith 
48015091d37SBarry Smith   Level: advanced
48115091d37SBarry Smith 
4821cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
48382bf6240SBarry Smith @*/
484d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
485d71ae5a4SJacob Faibussowitsch {
48682bf6240SBarry Smith   PetscFunctionBegin;
4870700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
488cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49082bf6240SBarry Smith }
49182bf6240SBarry Smith 
492ac226902SBarry Smith /*@C
49382bf6240SBarry Smith   TSPseudoSetTimeStep - Sets the user-defined routine to be
49482bf6240SBarry Smith   called at each pseudo-timestep to update the timestep.
49582bf6240SBarry Smith 
496c3339decSBarry Smith   Logically Collective
49715091d37SBarry Smith 
49882bf6240SBarry Smith   Input Parameters:
49915091d37SBarry Smith + ts  - timestep context
50082bf6240SBarry Smith . dt  - function to compute timestep
5012fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`)
50282bf6240SBarry Smith 
50320f4b53cSBarry Smith   Calling sequence of `dt`:
50414d0ab18SJacob Faibussowitsch + ts    - the `TS` context
50514d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep
50614d0ab18SJacob Faibussowitsch - ctx   - [optional] user-defined context
50782bf6240SBarry Smith 
508bcf0153eSBarry Smith   Level: intermediate
50982bf6240SBarry Smith 
510bcf0153eSBarry Smith   Notes:
511bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoComputeTimeStep()`
512bcf0153eSBarry Smith   during the timestepping process.
513bcf0153eSBarry Smith 
514bcf0153eSBarry Smith   If not set then `TSPseudoTimeStepDefault()` is automatically used
515bcf0153eSBarry Smith 
5161cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
51782bf6240SBarry Smith @*/
51814d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx)
519d71ae5a4SJacob Faibussowitsch {
52082bf6240SBarry Smith   PetscFunctionBegin;
5210700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
522cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52482bf6240SBarry Smith }
52582bf6240SBarry Smith 
52682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
52782bf6240SBarry Smith 
528ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/
529d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx)
530d71ae5a4SJacob Faibussowitsch {
531be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
53282bf6240SBarry Smith 
53382bf6240SBarry Smith   PetscFunctionBegin;
53482bf6240SBarry Smith   pseudo->verify    = dt;
53582bf6240SBarry Smith   pseudo->verifyctx = ctx;
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53782bf6240SBarry Smith }
53882bf6240SBarry Smith 
539d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc)
540d71ae5a4SJacob Faibussowitsch {
5414bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
54282bf6240SBarry Smith 
54382bf6240SBarry Smith   PetscFunctionBegin;
54482bf6240SBarry Smith   pseudo->dt_increment = inc;
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54682bf6240SBarry Smith }
54782bf6240SBarry Smith 
548d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt)
549d71ae5a4SJacob Faibussowitsch {
55086534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
55186534af1SJed Brown 
55286534af1SJed Brown   PetscFunctionBegin;
55386534af1SJed Brown   pseudo->dt_max = maxdt;
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55586534af1SJed Brown }
55686534af1SJed Brown 
557d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
558d71ae5a4SJacob Faibussowitsch {
5594bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
56082bf6240SBarry Smith 
56182bf6240SBarry Smith   PetscFunctionBegin;
5624bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56482bf6240SBarry Smith }
56582bf6240SBarry Smith 
5666849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/
567d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx)
568d71ae5a4SJacob Faibussowitsch {
5694bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
57082bf6240SBarry Smith 
57182bf6240SBarry Smith   PetscFunctionBegin;
57282bf6240SBarry Smith   pseudo->dt    = dt;
57382bf6240SBarry Smith   pseudo->dtctx = ctx;
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57582bf6240SBarry Smith }
57682bf6240SBarry Smith 
57710e6a065SJed Brown /*MC
5781d27aa22SBarry Smith   TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97`
57982bf6240SBarry Smith 
58010e6a065SJed Brown   This method solves equations of the form
58110e6a065SJed Brown 
5821d27aa22SBarry Smith   $$
5831d27aa22SBarry Smith   F(X,Xdot) = 0
5841d27aa22SBarry Smith   $$
58510e6a065SJed Brown 
58610e6a065SJed Brown   for steady state using the iteration
58710e6a065SJed Brown 
5881d27aa22SBarry Smith   $$
58937fdd005SBarry Smith   [G'] S = -F(X,0)
59037fdd005SBarry Smith   X += S
5911d27aa22SBarry Smith   $$
59210e6a065SJed Brown 
59310e6a065SJed Brown   where
59410e6a065SJed Brown 
5951d27aa22SBarry Smith   $$
5961d27aa22SBarry Smith   G(Y) = F(Y,(Y-X)/dt)
5971d27aa22SBarry Smith   $$
59810e6a065SJed Brown 
5996f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6006f2d6a7bSJed Brown   state".  See note below.
60110e6a065SJed Brown 
60293a54799SPierre Jolivet   In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`.
603fb59f417SJames Wright   It determines the next timestep via
604fb59f417SJames Wright 
605fb59f417SJames Wright   $$
606fb59f417SJames Wright   dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert}
607fb59f417SJames Wright   $$
608fb59f417SJames Wright 
609fb59f417SJames Wright   where $r$ is an additional growth factor (set by `-ts_pseudo_increment`).
610fb59f417SJames Wright   An alternative formulation is also available that uses the initial timestep and function norm.
611fb59f417SJames Wright 
612fb59f417SJames Wright   $$
613fb59f417SJames Wright   dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert}
614fb59f417SJames Wright   $$
615fb59f417SJames Wright 
616fb59f417SJames Wright   This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`.
617fb59f417SJames Wright   For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`.
618fb59f417SJames Wright 
619bcf0153eSBarry Smith   Options Database Keys:
62010e6a065SJed Brown +  -ts_pseudo_increment <real>                     - ratio of increase dt
6213118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
622fb59f417SJames Wright .  -ts_pseudo_max_dt                               - Maximum dt for adaptive timestepping algorithm
623fb59f417SJames Wright .  -ts_pseudo_monitor                              - Monitor convergence of the function norm
624fb59f417SJames Wright .  -ts_pseudo_fatol <atol>                         - stop iterating when the function norm is less than `atol`
625fb59f417SJames Wright -  -ts_pseudo_frtol <rtol>                         - stop iterating when the function norm divided by the initial function norm is less than `rtol`
62610e6a065SJed Brown 
62710e6a065SJed Brown   Level: beginner
62810e6a065SJed Brown 
62910e6a065SJed Brown   Notes:
6306f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6316f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6326f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6336f2d6a7bSJed Brown 
6341d27aa22SBarry Smith   $$
6351d27aa22SBarry Smith   Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6361d27aa22SBarry Smith   $$
6376f2d6a7bSJed Brown 
638*eb636060SBarry Smith   The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it
639*eb636060SBarry Smith   is $ dF/dX + shift*1/dt $. Hence  still contains the $ 1/dt $ term so the Jacobian is not simply the
640*eb636060SBarry Smith   Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$.
641*eb636060SBarry Smith 
6426f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6436f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6446f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
645fb59f417SJames Wright   By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers.
646fb59f417SJames Wright   Setting the `SNESType` via `-snes_type` will override this default setting.
64710e6a065SJed Brown 
6481cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`
64910e6a065SJed Brown M*/
650d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
651d71ae5a4SJacob Faibussowitsch {
6527bf11e45SBarry Smith   TS_Pseudo *pseudo;
653193ac0bcSJed Brown   SNES       snes;
65419fd82e9SBarry Smith   SNESType   stype;
6552d3f70b5SBarry Smith 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
657277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
658000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
659000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
660000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
661000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
662000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
6630f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
6640f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
6652ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
666825ab935SBarry Smith   ts->usessnes            = PETSC_TRUE;
6677bf11e45SBarry Smith 
6689566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
6699566063dSJacob Faibussowitsch   PetscCall(SNESGetType(snes, &stype));
6709566063dSJacob Faibussowitsch   if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY));
6712d3f70b5SBarry Smith 
6724dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pseudo));
6737bf11e45SBarry Smith   ts->data = (void *)pseudo;
6742d3f70b5SBarry Smith 
675be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
676be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
67728aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6784bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
679193ac0bcSJed Brown   pseudo->fnorm                        = -1;
680be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
681be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
6823118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
6833118ae5eSBarry Smith   pseudo->fatol = 1.e-25;
6843118ae5eSBarry Smith   pseudo->frtol = 1.e-5;
6853118ae5eSBarry Smith #else
6863118ae5eSBarry Smith   pseudo->fatol = 1.e-50;
6873118ae5eSBarry Smith   pseudo->frtol = 1.e-12;
6883118ae5eSBarry Smith #endif
6899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo));
6909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo));
6919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo));
6929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo));
6939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo));
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6952d3f70b5SBarry Smith }
6962d3f70b5SBarry Smith 
69782bf6240SBarry Smith /*@C
698bcf0153eSBarry Smith   TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.  Use with `TSPseudoSetTimeStep()`.
69928aa8177SBarry Smith 
700cc4c1da9SBarry Smith   Collective, No Fortran Support
70115091d37SBarry Smith 
70228aa8177SBarry Smith   Input Parameters:
703a2b725a8SWilliam Gropp + ts    - the timestep context
704a2b725a8SWilliam Gropp - dtctx - unused timestep context
70528aa8177SBarry Smith 
70682bf6240SBarry Smith   Output Parameter:
70782bf6240SBarry Smith . newdt - the timestep to use for the next step
70828aa8177SBarry Smith 
70915091d37SBarry Smith   Level: advanced
71015091d37SBarry Smith 
7111cc06b55SBarry Smith .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO`
71228aa8177SBarry Smith @*/
713d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx)
714d71ae5a4SJacob Faibussowitsch {
71582bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
716be5899b3SLisandro Dalcin   PetscReal  inc    = pseudo->dt_increment;
71728aa8177SBarry Smith 
7183a40ed3dSBarry Smith   PetscFunctionBegin;
71998473560SBarry Smith   if (pseudo->fnorm < 0.0) {
7209566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
7219566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
722*eb636060SBarry Smith     PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate));
7239566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
72498473560SBarry Smith   }
725be5899b3SLisandro Dalcin   if (pseudo->fnorm_initial < 0) {
72682bf6240SBarry Smith     /* first time through so compute initial function norm */
727cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial  = pseudo->fnorm;
728be5899b3SLisandro Dalcin     pseudo->fnorm_previous = pseudo->fnorm;
72982bf6240SBarry Smith   }
730bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step;
731bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm;
732be5899b3SLisandro Dalcin   else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm;
73386534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max);
73482bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73628aa8177SBarry Smith }
737