xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 940ebdd4d53c9bd2466ad09f7a10e4b0eba2bd47)
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 
6c558785fSJames Wright #define TSADAPTTSPSEUDO "tspseudo"
7c558785fSJames Wright 
82d3f70b5SBarry Smith typedef struct {
92d3f70b5SBarry Smith   Vec update; /* work vector where new solution is formed */
102d3f70b5SBarry Smith   Vec func;   /* work vector where F(t[i],u[i]) is stored */
116f2d6a7bSJed Brown   Vec xdot;   /* work vector for time derivative of state */
122d3f70b5SBarry Smith 
132d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
142d3f70b5SBarry Smith 
156849ba73SBarry Smith   PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */
162d3f70b5SBarry Smith   void *dtctx;
17ace3abfcSBarry Smith   PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */
187bf11e45SBarry Smith   void *verifyctx;
192d3f70b5SBarry Smith 
20cdbf8f93SLisandro Dalcin   PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */
2187828ca2SBarry Smith   PetscReal fnorm_previous;
2228aa8177SBarry Smith 
23cdbf8f93SLisandro Dalcin   PetscReal dt_initial;   /* initial time-step */
2487828ca2SBarry Smith   PetscReal dt_increment; /* scaling that dt is incremented each time-step */
2586534af1SJed Brown   PetscReal dt_max;       /* maximum time step */
26ace3abfcSBarry Smith   PetscBool increment_dt_from_initial_dt;
273118ae5eSBarry Smith   PetscReal fatol, frtol;
28eb636060SBarry Smith 
29eb636060SBarry Smith   PetscObjectState Xstate; /* state of input vector to TSComputeIFunction() with 0 Xdot */
30*940ebdd4SJames Wright 
31*940ebdd4SJames Wright   TSStepStatus status;
327bf11e45SBarry Smith } TS_Pseudo;
332d3f70b5SBarry Smith 
340fe2f5c2SJames Wright /*@C
350fe2f5c2SJames Wright   TSPseudoComputeFunction - Compute nonlinear residual for pseudo-timestepping
360fe2f5c2SJames Wright 
370fe2f5c2SJames Wright   This computes the residual for $\dot U = 0$, i.e. $F(U, 0)$ for the IFunction.
380fe2f5c2SJames Wright 
390fe2f5c2SJames Wright   Collective
400fe2f5c2SJames Wright 
410fe2f5c2SJames Wright   Input Parameters:
420fe2f5c2SJames Wright + ts       - the timestep context
430fe2f5c2SJames Wright - solution - the solution vector
440fe2f5c2SJames Wright 
450fe2f5c2SJames Wright   Output Parameter:
460fe2f5c2SJames Wright + residual - the nonlinear residual
470fe2f5c2SJames Wright - fnorm    - the norm of the nonlinear residual
480fe2f5c2SJames Wright 
490fe2f5c2SJames Wright   Level: advanced
500fe2f5c2SJames Wright 
510fe2f5c2SJames Wright   Note:
520fe2f5c2SJames Wright   `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.
530fe2f5c2SJames Wright 
540fe2f5c2SJames Wright   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.
550fe2f5c2SJames Wright 
560fe2f5c2SJames Wright .seealso: [](ch_ts), `TSPSEUDO`
570fe2f5c2SJames Wright @*/
580fe2f5c2SJames Wright PetscErrorCode TSPseudoComputeFunction(TS ts, Vec solution, Vec *residual, PetscReal *fnorm)
590fe2f5c2SJames Wright {
600fe2f5c2SJames Wright   TS_Pseudo       *pseudo = (TS_Pseudo *)ts->data;
610fe2f5c2SJames Wright   PetscObjectState Xstate;
620fe2f5c2SJames Wright 
630fe2f5c2SJames Wright   PetscFunctionBegin;
640fe2f5c2SJames Wright   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
650fe2f5c2SJames Wright   PetscValidHeaderSpecific(solution, VEC_CLASSID, 2);
660fe2f5c2SJames Wright   if (residual) PetscValidHeaderSpecific(*residual, VEC_CLASSID, 3);
670fe2f5c2SJames Wright   if (fnorm) PetscAssertPointer(fnorm, 4);
680fe2f5c2SJames Wright 
690fe2f5c2SJames Wright   PetscCall(PetscObjectStateGet((PetscObject)solution, &Xstate));
700fe2f5c2SJames Wright   if (Xstate != pseudo->Xstate || pseudo->fnorm < 0) {
710fe2f5c2SJames Wright     PetscCall(VecZeroEntries(pseudo->xdot));
720fe2f5c2SJames Wright     PetscCall(TSComputeIFunction(ts, ts->ptime, solution, pseudo->xdot, pseudo->func, PETSC_FALSE));
730fe2f5c2SJames Wright     pseudo->Xstate = Xstate;
740fe2f5c2SJames Wright     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
750fe2f5c2SJames Wright   }
760fe2f5c2SJames Wright   if (residual) *residual = pseudo->func;
770fe2f5c2SJames Wright   if (fnorm) *fnorm = pseudo->fnorm;
780fe2f5c2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
790fe2f5c2SJames Wright }
800fe2f5c2SJames Wright 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts)
82d71ae5a4SJacob Faibussowitsch {
83277b19d0SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
84*940ebdd4SJames Wright   PetscInt   nits, lits, rejections = 0;
85*940ebdd4SJames Wright   PetscBool  accept;
86be5899b3SLisandro Dalcin   PetscReal  next_time_step = ts->time_step;
87c558785fSJames Wright   TSAdapt    adapt;
882d3f70b5SBarry Smith 
893a40ed3dSBarry Smith   PetscFunctionBegin;
90eb636060SBarry Smith   if (ts->steps == 0) {
91eb636060SBarry Smith     pseudo->dt_initial = ts->time_step;
92eb636060SBarry Smith     PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* in all future updates pseudo->update already contains the current time solution */
93eb636060SBarry Smith   }
94*940ebdd4SJames Wright 
95*940ebdd4SJames Wright   pseudo->status = TS_STEP_INCOMPLETE;
96*940ebdd4SJames Wright   while (!ts->reason && pseudo->status != TS_STEP_COMPLETE) {
97*940ebdd4SJames Wright     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 */
989566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime + ts->time_step));
999566063dSJacob Faibussowitsch     PetscCall(SNESSolve(ts->snes, NULL, pseudo->update));
1009566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1019566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1029371c9d4SSatish Balay     ts->snes_its += nits;
1039371c9d4SSatish Balay     ts->ksp_its += lits;
104c558785fSJames Wright     pseudo->fnorm = -1; /* The current norm is no longer valid */
105f4f49eeaSPierre Jolivet     PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &pseudo->update));
106c558785fSJames Wright     PetscCall(TSGetAdapt(ts, &adapt));
107*940ebdd4SJames Wright     PetscCall(TSAdaptCheckStage(adapt, ts, ts->ptime + ts->time_step, pseudo->update, &accept));
108*940ebdd4SJames Wright     if (!accept) goto reject_step;
109be5899b3SLisandro Dalcin 
110*940ebdd4SJames Wright     pseudo->status = TS_STEP_PENDING;
111*940ebdd4SJames Wright     PetscCall(VecCopy(pseudo->update, ts->vec_sol));
112*940ebdd4SJames Wright     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
113*940ebdd4SJames Wright     pseudo->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
114*940ebdd4SJames Wright     if (!accept) {
115*940ebdd4SJames Wright       ts->time_step = next_time_step;
116*940ebdd4SJames Wright       goto reject_step;
117*940ebdd4SJames Wright     }
118be5899b3SLisandro Dalcin     ts->ptime += ts->time_step;
119be5899b3SLisandro Dalcin     ts->time_step = next_time_step;
120*940ebdd4SJames Wright     break;
121be5899b3SLisandro Dalcin 
122*940ebdd4SJames Wright   reject_step:
123*940ebdd4SJames Wright     ts->reject++;
124*940ebdd4SJames Wright     accept = PETSC_FALSE;
125*940ebdd4SJames Wright     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
126*940ebdd4SJames Wright       ts->reason = TS_DIVERGED_STEP_REJECTED;
127*940ebdd4SJames Wright       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
128*940ebdd4SJames Wright       PetscFunctionReturn(PETSC_SUCCESS);
129*940ebdd4SJames Wright     }
130*940ebdd4SJames Wright   }
131*940ebdd4SJames Wright 
132*940ebdd4SJames Wright   // Check solution convergence
1330fe2f5c2SJames Wright   PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, NULL));
134eb636060SBarry Smith 
1353118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1363118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
13763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol));
1383ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1393118ae5eSBarry Smith   }
1403118ae5eSBarry Smith   if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) {
1413118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
14263a3b9bcSJacob 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));
1433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1443118ae5eSBarry Smith   }
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1462d3f70b5SBarry Smith }
1472d3f70b5SBarry Smith 
148d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts)
149d71ae5a4SJacob Faibussowitsch {
1507bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1512d3f70b5SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->update));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->func));
1559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->xdot));
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1572d3f70b5SBarry Smith }
1582d3f70b5SBarry Smith 
159d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts)
160d71ae5a4SJacob Faibussowitsch {
161277b19d0SLisandro Dalcin   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(TSReset_Pseudo(ts));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
1649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL));
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL));
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL));
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL));
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170277b19d0SLisandro Dalcin }
1712d3f70b5SBarry Smith 
172d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot)
173d71ae5a4SJacob Faibussowitsch {
1746f2d6a7bSJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1752d3f70b5SBarry Smith 
1763a40ed3dSBarry Smith   PetscFunctionBegin;
1776f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1792d3f70b5SBarry Smith }
1802d3f70b5SBarry Smith 
1816f2d6a7bSJed Brown /*
1826f2d6a7bSJed Brown     The transient residual is
1836f2d6a7bSJed Brown 
1846f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
1856f2d6a7bSJed Brown 
1866f2d6a7bSJed Brown     or for ODE,
1876f2d6a7bSJed Brown 
1886f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
1896f2d6a7bSJed Brown 
1906f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
1916f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
1926f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
1936f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
1946f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
195eb636060SBarry Smith     algorithm, it only takes this one full Newton step with the steady state
1966f2d6a7bSJed Brown     residual, and then advances to the next time step.
197eb636060SBarry Smith 
198eb636060SBarry Smith     This routine avoids a repeated identical call to TSComputeRHSFunction() when that result
199eb636060SBarry Smith     is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo()
200eb636060SBarry Smith 
2016f2d6a7bSJed Brown */
202d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts)
203d71ae5a4SJacob Faibussowitsch {
204eb636060SBarry Smith   TSIFunctionFn    *ifunction;
205eb636060SBarry Smith   TSRHSFunctionFn  *rhsfunction;
206eb636060SBarry Smith   void             *ctx;
207eb636060SBarry Smith   DM                dm;
208eb636060SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo *)ts->data;
20981775c15SStefano Zampini   const PetscScalar mdt    = 1.0 / ts->time_step;
210eb636060SBarry Smith   PetscBool         KSPSNES;
211eb636060SBarry Smith   PetscObjectState  Xstate;
2122d3f70b5SBarry Smith 
2133a40ed3dSBarry Smith   PetscFunctionBegin;
214eb636060SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
215eb636060SBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
216eb636060SBarry Smith   PetscValidHeaderSpecific(Y, VEC_CLASSID, 3);
217eb636060SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 4);
218eb636060SBarry Smith 
219eb636060SBarry Smith   PetscCall(TSGetDM(ts, &dm));
220eb636060SBarry Smith   PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx));
221eb636060SBarry Smith   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
222eb636060SBarry Smith   PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()");
223eb636060SBarry Smith 
224eb636060SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESKSPONLY, &KSPSNES));
225eb636060SBarry Smith   PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate));
226eb636060SBarry Smith   if (Xstate != pseudo->Xstate || ifunction || !KSPSNES) {
227c558785fSJames Wright     // X == ts->vec_sol for first SNES iteration, so pseudo->xdot == 0
22881775c15SStefano Zampini     PetscCall(VecAXPBYPCZ(pseudo->xdot, -mdt, mdt, 0, ts->vec_sol, X));
22981775c15SStefano Zampini     PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, pseudo->xdot, Y, PETSC_FALSE));
230eb636060SBarry Smith   } else {
231eb636060SBarry Smith     /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */
232eb636060SBarry Smith     /* note that pseudo->func contains the negation of TSComputeRHSFunction() */
23381775c15SStefano Zampini     PetscCall(VecWAXPY(Y, 1, pseudo->func, pseudo->xdot));
234eb636060SBarry Smith   }
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2366f2d6a7bSJed Brown }
2372d3f70b5SBarry Smith 
2386f2d6a7bSJed Brown /*
2396f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2406f2d6a7bSJed Brown 
2416f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2426f2d6a7bSJed Brown 
2436f2d6a7bSJed Brown     and for ODE:
2446f2d6a7bSJed Brown 
2456f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2466f2d6a7bSJed Brown */
247d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts)
248d71ae5a4SJacob Faibussowitsch {
2496f2d6a7bSJed Brown   Vec Xdot;
2506f2d6a7bSJed Brown 
2516f2d6a7bSJed Brown   PetscFunctionBegin;
25281775c15SStefano Zampini   /* Xdot has already been computed in SNESTSFormFunction_Pseudo (SNES guarantees this) */
2539566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
2549566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE));
2553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2562d3f70b5SBarry Smith }
2572d3f70b5SBarry Smith 
258d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts)
259d71ae5a4SJacob Faibussowitsch {
2607bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
2612d3f70b5SBarry Smith 
2623a40ed3dSBarry Smith   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update));
2649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func));
2659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2672d3f70b5SBarry Smith }
26881775c15SStefano Zampini 
269d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy)
270d71ae5a4SJacob Faibussowitsch {
2717bf11e45SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
272ce94432eSBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
2732d3f70b5SBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
2750fe2f5c2SJames Wright   PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, NULL));
2769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
27763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm));
2789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2802d3f70b5SBarry Smith }
2812d3f70b5SBarry Smith 
282ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject)
283d71ae5a4SJacob Faibussowitsch {
2844bbc92c1SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
285ace3abfcSBarry Smith   PetscBool   flg    = PETSC_FALSE;
286649052a6SBarry Smith   PetscViewer viewer;
2872d3f70b5SBarry Smith 
2883a40ed3dSBarry Smith   PetscFunctionBegin;
289d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options");
2909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL));
2912d3f70b5SBarry Smith   if (flg) {
2929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer));
29349abdd8aSBarry Smith     PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
29428aa8177SBarry Smith   }
295be5899b3SLisandro Dalcin   flg = pseudo->increment_dt_from_initial_dt;
2969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL));
297be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
2989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL));
2999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL));
3009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL));
3019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL));
302d0609cedSBarry Smith   PetscOptionsHeadEnd();
3033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3042d3f70b5SBarry Smith }
3052d3f70b5SBarry Smith 
306d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer)
307d71ae5a4SJacob Faibussowitsch {
3083118ae5eSBarry Smith   PetscBool isascii;
309d52bd9f3SBarry Smith 
3103a40ed3dSBarry Smith   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3123118ae5eSBarry Smith   if (isascii) {
3133118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
3149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  frtol - relative tolerance in function value %g\n", (double)pseudo->frtol));
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol));
3169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_initial - initial timestep %g\n", (double)pseudo->dt_initial));
3179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment));
3189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_max - maximum time %g\n", (double)pseudo->dt_max));
3193118ae5eSBarry Smith   }
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3212d3f70b5SBarry Smith }
3222d3f70b5SBarry Smith 
323ac226902SBarry Smith /*@C
324c558785fSJames Wright   TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.  Use with `TSPseudoSetTimeStep()`.
325c558785fSJames Wright 
326c558785fSJames Wright   Collective, No Fortran Support
327c558785fSJames Wright 
328c558785fSJames Wright   Input Parameters:
329c558785fSJames Wright + ts    - the timestep context
330c558785fSJames Wright - dtctx - unused timestep context
331c558785fSJames Wright 
332c558785fSJames Wright   Output Parameter:
333c558785fSJames Wright . newdt - the timestep to use for the next step
334c558785fSJames Wright 
335c558785fSJames Wright   Level: advanced
336c558785fSJames Wright 
337c558785fSJames Wright .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO`
338c558785fSJames Wright @*/
339c558785fSJames Wright PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx)
340c558785fSJames Wright {
341c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
3420fe2f5c2SJames Wright   PetscReal  inc    = pseudo->dt_increment, fnorm;
343c558785fSJames Wright 
344c558785fSJames Wright   PetscFunctionBegin;
3450fe2f5c2SJames Wright   PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, &fnorm));
346c558785fSJames Wright   if (pseudo->fnorm_initial < 0) {
347c558785fSJames Wright     /* first time through so compute initial function norm */
3480fe2f5c2SJames Wright     pseudo->fnorm_initial  = fnorm;
3490fe2f5c2SJames Wright     pseudo->fnorm_previous = fnorm;
350c558785fSJames Wright   }
3510fe2f5c2SJames Wright   if (fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step;
3520fe2f5c2SJames Wright   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / fnorm;
3530fe2f5c2SJames Wright   else *newdt = inc * ts->time_step * pseudo->fnorm_previous / fnorm;
354c558785fSJames Wright   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max);
3550fe2f5c2SJames Wright   pseudo->fnorm_previous = fnorm;
356c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
357c558785fSJames Wright }
358c558785fSJames Wright 
359c558785fSJames Wright static PetscErrorCode TSAdaptChoose_TSPseudo(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
360c558785fSJames Wright {
361c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
362c558785fSJames Wright 
363c558785fSJames Wright   PetscFunctionBegin;
364c558785fSJames Wright   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
365c558785fSJames Wright   PetscCall((*pseudo->dt)(ts, next_h, pseudo->dtctx));
366c558785fSJames Wright   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
367c558785fSJames Wright 
368c558785fSJames Wright   *next_sc = 0;
369c558785fSJames Wright   *wlte    = -1; /* Weighted local truncation error was not evaluated */
370c558785fSJames Wright   *wltea   = -1; /* Weighted absolute local truncation error was not evaluated */
371c558785fSJames Wright   *wlter   = -1; /* Weighted relative local truncation error was not evaluated */
372c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
373c558785fSJames Wright }
374c558785fSJames Wright 
375c558785fSJames Wright static PetscErrorCode TSAdaptCheckStage_TSPseudo(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
376c558785fSJames Wright {
377c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
378c558785fSJames Wright 
379c558785fSJames Wright   PetscFunctionBegin;
380c558785fSJames Wright   if (pseudo->verify) {
381c558785fSJames Wright     PetscReal dt;
382c558785fSJames Wright     PetscCall(TSGetTimeStep(ts, &dt));
383c558785fSJames Wright     PetscCall((*pseudo->verify)(ts, Y, pseudo->verifyctx, &dt, accept));
384c558785fSJames Wright     PetscCall(TSSetTimeStep(ts, dt));
385c558785fSJames Wright   }
386c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
387c558785fSJames Wright }
388c558785fSJames Wright 
389c558785fSJames Wright /*MC
390c558785fSJames Wright   TSADAPTTSPSEUDO - TSPseudo adaptive controller for time stepping
391c558785fSJames Wright 
392c558785fSJames Wright   Level: developer
393c558785fSJames Wright 
394c558785fSJames Wright   Note:
395c558785fSJames Wright   This is only meant to be used with `TSPSEUDO` time integrator.
396c558785fSJames Wright 
397c558785fSJames Wright .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`, `TSPSEUDO`
398c558785fSJames Wright M*/
399c558785fSJames Wright static PetscErrorCode TSAdaptCreate_TSPseudo(TSAdapt adapt)
400c558785fSJames Wright {
401c558785fSJames Wright   PetscFunctionBegin;
402c558785fSJames Wright   adapt->ops->choose = TSAdaptChoose_TSPseudo;
403c558785fSJames Wright   adapt->checkstage  = TSAdaptCheckStage_TSPseudo;
404c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
405c558785fSJames Wright }
406c558785fSJames Wright 
407c558785fSJames Wright /*@
408c558785fSJames Wright   TSPseudoComputeTimeStep - Computes the next timestep for a currently running
409c558785fSJames Wright   pseudo-timestepping process.
410c558785fSJames Wright 
411c558785fSJames Wright   Collective
412c558785fSJames Wright 
413c558785fSJames Wright   Input Parameter:
414c558785fSJames Wright . ts - timestep context
415c558785fSJames Wright 
416c558785fSJames Wright   Output Parameter:
417c558785fSJames Wright . dt - newly computed timestep
418c558785fSJames Wright 
419c558785fSJames Wright   Level: developer
420c558785fSJames Wright 
421c558785fSJames Wright   Note:
422c558785fSJames Wright   The routine to be called here to compute the timestep should be
423c558785fSJames Wright   set by calling `TSPseudoSetTimeStep()`.
424c558785fSJames Wright 
425c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
426c558785fSJames Wright @*/
427c558785fSJames Wright PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt)
428c558785fSJames Wright {
429c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
430c558785fSJames Wright 
431c558785fSJames Wright   PetscFunctionBegin;
432c558785fSJames Wright   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
433c558785fSJames Wright   PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx));
434c558785fSJames Wright   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
435c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
436c558785fSJames Wright }
437c558785fSJames Wright 
438c558785fSJames Wright /*@C
439c558785fSJames Wright   TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
440c558785fSJames Wright 
441c558785fSJames Wright   Collective, No Fortran Support
442c558785fSJames Wright 
443c558785fSJames Wright   Input Parameters:
444c558785fSJames Wright + ts     - the timestep context
445c558785fSJames Wright . dtctx  - unused timestep context
446c558785fSJames Wright - update - latest solution vector
447c558785fSJames Wright 
448c558785fSJames Wright   Output Parameters:
449c558785fSJames Wright + newdt - the timestep to use for the next step
450c558785fSJames Wright - flag  - flag indicating whether the last time step was acceptable
451c558785fSJames Wright 
452c558785fSJames Wright   Level: advanced
453c558785fSJames Wright 
454c558785fSJames Wright   Note:
455c558785fSJames Wright   This routine always returns a flag of 1, indicating an acceptable
456c558785fSJames Wright   timestep.
457c558785fSJames Wright 
458c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
459c558785fSJames Wright @*/
460c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag)
461c558785fSJames Wright {
462c558785fSJames Wright   PetscFunctionBegin;
463c558785fSJames Wright   // NOTE: This function was never used
464c558785fSJames Wright   *flag = PETSC_TRUE;
465c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
466c558785fSJames Wright }
467c558785fSJames Wright 
468c558785fSJames Wright /*@
469c558785fSJames Wright   TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
470c558785fSJames Wright 
471c558785fSJames Wright   Collective
472c558785fSJames Wright 
473c558785fSJames Wright   Input Parameters:
474c558785fSJames Wright + ts     - timestep context
475c558785fSJames Wright - update - latest solution vector
476c558785fSJames Wright 
477c558785fSJames Wright   Output Parameters:
478c558785fSJames Wright + dt   - newly computed timestep (if it had to shrink)
479c558785fSJames Wright - flag - indicates if current timestep was ok
480c558785fSJames Wright 
481c558785fSJames Wright   Level: advanced
482c558785fSJames Wright 
483c558785fSJames Wright   Notes:
484c558785fSJames Wright   The routine to be called here to compute the timestep should be
485c558785fSJames Wright   set by calling `TSPseudoSetVerifyTimeStep()`.
486c558785fSJames Wright 
487c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
488c558785fSJames Wright @*/
489c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag)
490c558785fSJames Wright {
491c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
492c558785fSJames Wright 
493c558785fSJames Wright   PetscFunctionBegin;
494c558785fSJames Wright   // NOTE: This function is never used
495c558785fSJames Wright   *flag = PETSC_TRUE;
496c558785fSJames Wright   if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag));
497c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
498c558785fSJames Wright }
499c558785fSJames Wright 
500c558785fSJames Wright /*@C
50182bf6240SBarry Smith   TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
50282bf6240SBarry Smith   last timestep.
50382bf6240SBarry Smith 
504c3339decSBarry Smith   Logically Collective
50515091d37SBarry Smith 
50682bf6240SBarry Smith   Input Parameters:
50715091d37SBarry Smith + ts  - timestep context
50882bf6240SBarry Smith . dt  - user-defined function to verify timestep
5092fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`)
51082bf6240SBarry Smith 
51120f4b53cSBarry Smith   Calling sequence of `func`:
5122fe279fdSBarry Smith + ts     - the time-step context
5132fe279fdSBarry Smith . update - latest solution vector
51414d0ab18SJacob Faibussowitsch . ctx    - [optional] user-defined timestep context
51582bf6240SBarry Smith . newdt  - the timestep to use for the next step
516a2b725a8SWilliam Gropp - flag   - flag indicating whether the last time step was acceptable
51782bf6240SBarry Smith 
518bcf0153eSBarry Smith   Level: advanced
519bcf0153eSBarry Smith 
520bcf0153eSBarry Smith   Note:
521bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoVerifyTimeStep()`
52282bf6240SBarry Smith   during the timestepping process.
52382bf6240SBarry Smith 
5241cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
52582bf6240SBarry Smith @*/
52614d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx)
527d71ae5a4SJacob Faibussowitsch {
52882bf6240SBarry Smith   PetscFunctionBegin;
5290700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
530cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx));
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53282bf6240SBarry Smith }
53382bf6240SBarry Smith 
53482bf6240SBarry Smith /*@
53582bf6240SBarry Smith   TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
5368d359177SBarry Smith   dt when using the TSPseudoTimeStepDefault() routine.
53782bf6240SBarry Smith 
538c3339decSBarry Smith   Logically Collective
539fee21e36SBarry Smith 
54015091d37SBarry Smith   Input Parameters:
54115091d37SBarry Smith + ts  - the timestep context
54215091d37SBarry Smith - inc - the scaling factor >= 1.0
54315091d37SBarry Smith 
54482bf6240SBarry Smith   Options Database Key:
54567b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment
54682bf6240SBarry Smith 
54715091d37SBarry Smith   Level: advanced
54815091d37SBarry Smith 
5491cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
55082bf6240SBarry Smith @*/
551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc)
552d71ae5a4SJacob Faibussowitsch {
55382bf6240SBarry Smith   PetscFunctionBegin;
5540700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
555c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts, inc, 2);
556cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc));
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55882bf6240SBarry Smith }
55982bf6240SBarry Smith 
56086534af1SJed Brown /*@
56186534af1SJed Brown   TSPseudoSetMaxTimeStep - Sets the maximum time step
5628d359177SBarry Smith   when using the TSPseudoTimeStepDefault() routine.
56386534af1SJed Brown 
564c3339decSBarry Smith   Logically Collective
56586534af1SJed Brown 
56686534af1SJed Brown   Input Parameters:
56786534af1SJed Brown + ts    - the timestep context
56886534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate
56986534af1SJed Brown 
57086534af1SJed Brown   Options Database Key:
57167b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt
57286534af1SJed Brown 
57386534af1SJed Brown   Level: advanced
57486534af1SJed Brown 
5751cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
57686534af1SJed Brown @*/
577d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt)
578d71ae5a4SJacob Faibussowitsch {
57986534af1SJed Brown   PetscFunctionBegin;
58086534af1SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
58186534af1SJed Brown   PetscValidLogicalCollectiveReal(ts, maxdt, 2);
582cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
5833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58486534af1SJed Brown }
58586534af1SJed Brown 
58682bf6240SBarry Smith /*@
58782bf6240SBarry Smith   TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
588b44f4de4SBarry 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.$
58982bf6240SBarry Smith 
590c3339decSBarry Smith   Logically Collective
59115091d37SBarry Smith 
59282bf6240SBarry Smith   Input Parameter:
59382bf6240SBarry Smith . ts - the timestep context
59482bf6240SBarry Smith 
59582bf6240SBarry Smith   Options Database Key:
596b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment
59782bf6240SBarry Smith 
59815091d37SBarry Smith   Level: advanced
59915091d37SBarry Smith 
6001cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
60182bf6240SBarry Smith @*/
602d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
603d71ae5a4SJacob Faibussowitsch {
60482bf6240SBarry Smith   PetscFunctionBegin;
6050700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
606cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts));
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60882bf6240SBarry Smith }
60982bf6240SBarry Smith 
610ac226902SBarry Smith /*@C
61182bf6240SBarry Smith   TSPseudoSetTimeStep - Sets the user-defined routine to be
61282bf6240SBarry Smith   called at each pseudo-timestep to update the timestep.
61382bf6240SBarry Smith 
614c3339decSBarry Smith   Logically Collective
61515091d37SBarry Smith 
61682bf6240SBarry Smith   Input Parameters:
61715091d37SBarry Smith + ts  - timestep context
61882bf6240SBarry Smith . dt  - function to compute timestep
6192fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`)
62082bf6240SBarry Smith 
62120f4b53cSBarry Smith   Calling sequence of `dt`:
62214d0ab18SJacob Faibussowitsch + ts    - the `TS` context
62314d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep
62414d0ab18SJacob Faibussowitsch - ctx   - [optional] user-defined context
62582bf6240SBarry Smith 
626bcf0153eSBarry Smith   Level: intermediate
62782bf6240SBarry Smith 
628bcf0153eSBarry Smith   Notes:
629bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoComputeTimeStep()`
630bcf0153eSBarry Smith   during the timestepping process.
631bcf0153eSBarry Smith 
632bcf0153eSBarry Smith   If not set then `TSPseudoTimeStepDefault()` is automatically used
633bcf0153eSBarry Smith 
6341cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
63582bf6240SBarry Smith @*/
63614d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx)
637d71ae5a4SJacob Faibussowitsch {
63882bf6240SBarry Smith   PetscFunctionBegin;
6390700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
640cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx));
6413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64282bf6240SBarry Smith }
64382bf6240SBarry Smith 
644ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/
645d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx)
646d71ae5a4SJacob Faibussowitsch {
647be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
64882bf6240SBarry Smith 
64982bf6240SBarry Smith   PetscFunctionBegin;
65082bf6240SBarry Smith   pseudo->verify    = dt;
65182bf6240SBarry Smith   pseudo->verifyctx = ctx;
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65382bf6240SBarry Smith }
65482bf6240SBarry Smith 
655d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc)
656d71ae5a4SJacob Faibussowitsch {
6574bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
65882bf6240SBarry Smith 
65982bf6240SBarry Smith   PetscFunctionBegin;
66082bf6240SBarry Smith   pseudo->dt_increment = inc;
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66282bf6240SBarry Smith }
66382bf6240SBarry Smith 
664d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt)
665d71ae5a4SJacob Faibussowitsch {
66686534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
66786534af1SJed Brown 
66886534af1SJed Brown   PetscFunctionBegin;
66986534af1SJed Brown   pseudo->dt_max = maxdt;
6703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67186534af1SJed Brown }
67286534af1SJed Brown 
673d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
674d71ae5a4SJacob Faibussowitsch {
6754bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
67682bf6240SBarry Smith 
67782bf6240SBarry Smith   PetscFunctionBegin;
6784bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68082bf6240SBarry Smith }
68182bf6240SBarry Smith 
6826849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/
683d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx)
684d71ae5a4SJacob Faibussowitsch {
6854bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
68682bf6240SBarry Smith 
68782bf6240SBarry Smith   PetscFunctionBegin;
68882bf6240SBarry Smith   pseudo->dt    = dt;
68982bf6240SBarry Smith   pseudo->dtctx = ctx;
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69182bf6240SBarry Smith }
69282bf6240SBarry Smith 
69310e6a065SJed Brown /*MC
6941d27aa22SBarry Smith   TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97`
69582bf6240SBarry Smith 
69610e6a065SJed Brown   This method solves equations of the form
69710e6a065SJed Brown 
6981d27aa22SBarry Smith   $$
6991d27aa22SBarry Smith   F(X,Xdot) = 0
7001d27aa22SBarry Smith   $$
70110e6a065SJed Brown 
70210e6a065SJed Brown   for steady state using the iteration
70310e6a065SJed Brown 
7041d27aa22SBarry Smith   $$
70537fdd005SBarry Smith   [G'] S = -F(X,0)
70637fdd005SBarry Smith   X += S
7071d27aa22SBarry Smith   $$
70810e6a065SJed Brown 
70910e6a065SJed Brown   where
71010e6a065SJed Brown 
7111d27aa22SBarry Smith   $$
7121d27aa22SBarry Smith   G(Y) = F(Y,(Y-X)/dt)
7131d27aa22SBarry Smith   $$
71410e6a065SJed Brown 
7156f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
7166f2d6a7bSJed Brown   state".  See note below.
71710e6a065SJed Brown 
71893a54799SPierre Jolivet   In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`.
719fb59f417SJames Wright   It determines the next timestep via
720fb59f417SJames Wright 
721fb59f417SJames Wright   $$
722fb59f417SJames Wright   dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert}
723fb59f417SJames Wright   $$
724fb59f417SJames Wright 
725fb59f417SJames Wright   where $r$ is an additional growth factor (set by `-ts_pseudo_increment`).
726fb59f417SJames Wright   An alternative formulation is also available that uses the initial timestep and function norm.
727fb59f417SJames Wright 
728fb59f417SJames Wright   $$
729fb59f417SJames Wright   dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert}
730fb59f417SJames Wright   $$
731fb59f417SJames Wright 
732fb59f417SJames Wright   This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`.
733fb59f417SJames Wright   For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`.
734fb59f417SJames Wright 
735bcf0153eSBarry Smith   Options Database Keys:
73610e6a065SJed Brown +  -ts_pseudo_increment <real>                     - ratio of increase dt
7373118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
738fb59f417SJames Wright .  -ts_pseudo_max_dt                               - Maximum dt for adaptive timestepping algorithm
739fb59f417SJames Wright .  -ts_pseudo_monitor                              - Monitor convergence of the function norm
740fb59f417SJames Wright .  -ts_pseudo_fatol <atol>                         - stop iterating when the function norm is less than `atol`
741fb59f417SJames Wright -  -ts_pseudo_frtol <rtol>                         - stop iterating when the function norm divided by the initial function norm is less than `rtol`
74210e6a065SJed Brown 
74310e6a065SJed Brown   Level: beginner
74410e6a065SJed Brown 
74510e6a065SJed Brown   Notes:
7466f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
7476f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
7486f2d6a7bSJed Brown   last step, on the first Newton iteration we have
7496f2d6a7bSJed Brown 
7501d27aa22SBarry Smith   $$
7511d27aa22SBarry Smith   Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
7521d27aa22SBarry Smith   $$
7536f2d6a7bSJed Brown 
754eb636060SBarry Smith   The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it
755eb636060SBarry Smith   is $ dF/dX + shift*1/dt $. Hence  still contains the $ 1/dt $ term so the Jacobian is not simply the
756eb636060SBarry Smith   Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$.
757eb636060SBarry Smith 
7586f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
7596f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
7606f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
761fb59f417SJames Wright   By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers.
762fb59f417SJames Wright   Setting the `SNESType` via `-snes_type` will override this default setting.
76310e6a065SJed Brown 
7641cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`
76510e6a065SJed Brown M*/
766d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
767d71ae5a4SJacob Faibussowitsch {
7687bf11e45SBarry Smith   TS_Pseudo *pseudo;
769193ac0bcSJed Brown   SNES       snes;
77019fd82e9SBarry Smith   SNESType   stype;
7712d3f70b5SBarry Smith 
7723a40ed3dSBarry Smith   PetscFunctionBegin;
773277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
774000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
775000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
776000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
777000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
778000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
7790f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
7800f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
781c558785fSJames Wright   ts->default_adapt_type  = TSADAPTTSPSEUDO;
782825ab935SBarry Smith   ts->usessnes            = PETSC_TRUE;
7837bf11e45SBarry Smith 
784c558785fSJames Wright   PetscCall(TSAdaptRegister(TSADAPTTSPSEUDO, TSAdaptCreate_TSPseudo));
785c558785fSJames Wright 
7869566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
7879566063dSJacob Faibussowitsch   PetscCall(SNESGetType(snes, &stype));
7889566063dSJacob Faibussowitsch   if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY));
7892d3f70b5SBarry Smith 
7904dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pseudo));
7917bf11e45SBarry Smith   ts->data = (void *)pseudo;
7922d3f70b5SBarry Smith 
793be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
794be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
79528aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7964bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
797193ac0bcSJed Brown   pseudo->fnorm                        = -1;
798be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
799be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
8003118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
8013118ae5eSBarry Smith   pseudo->fatol = 1.e-25;
8023118ae5eSBarry Smith   pseudo->frtol = 1.e-5;
8033118ae5eSBarry Smith #else
8043118ae5eSBarry Smith   pseudo->fatol = 1.e-50;
8053118ae5eSBarry Smith   pseudo->frtol = 1.e-12;
8063118ae5eSBarry Smith #endif
8079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo));
8089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo));
8099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo));
8109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo));
8119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo));
8123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8132d3f70b5SBarry Smith }
814