xref: /petsc/src/ts/impls/pseudo/posindep.c (revision c558785f0dca1571fa7f63d87fda0eec410d3843)
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 
6*c558785fSJames Wright #define TSADAPTTSPSEUDO "tspseudo"
7*c558785fSJames 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 */
307bf11e45SBarry Smith } TS_Pseudo;
312d3f70b5SBarry Smith 
32d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts)
33d71ae5a4SJacob Faibussowitsch {
34277b19d0SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
35be5899b3SLisandro Dalcin   PetscInt   nits, lits, reject;
36cdbf8f93SLisandro Dalcin   PetscBool  stepok;
37be5899b3SLisandro Dalcin   PetscReal  next_time_step = ts->time_step;
38*c558785fSJames Wright   TSAdapt    adapt;
392d3f70b5SBarry Smith 
403a40ed3dSBarry Smith   PetscFunctionBegin;
41eb636060SBarry Smith   if (ts->steps == 0) {
42eb636060SBarry Smith     pseudo->dt_initial = ts->time_step;
43eb636060SBarry Smith     PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* in all future updates pseudo->update already contains the current time solution */
44eb636060SBarry Smith   }
45cdbf8f93SLisandro Dalcin   for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) {
46cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
47*c558785fSJames Wright     if (reject) PetscCall(VecCopy(ts->vec_sol0, pseudo->update)); /* need to copy the solution because we got a rejection and pseudo->update contains intermediate computation */
489566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime + ts->time_step));
499566063dSJacob Faibussowitsch     PetscCall(SNESSolve(ts->snes, NULL, pseudo->update));
509566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(ts->snes, &nits));
519566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
529371c9d4SSatish Balay     ts->snes_its += nits;
539371c9d4SSatish Balay     ts->ksp_its += lits;
54*c558785fSJames Wright     pseudo->fnorm = -1; /* The current norm is no longer valid */
55f4f49eeaSPierre Jolivet     PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &pseudo->update));
56*c558785fSJames Wright     PetscCall(TSGetAdapt(ts, &adapt));
57*c558785fSJames Wright     PetscCall(TSAdaptCheckStage(adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok));
589371c9d4SSatish Balay     if (!stepok) {
599371c9d4SSatish Balay       next_time_step = ts->time_step;
609371c9d4SSatish Balay       continue;
619371c9d4SSatish Balay     }
62*c558785fSJames Wright     PetscCall(VecCopy(pseudo->update, ts->vec_sol));
63*c558785fSJames Wright     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &stepok));
64cdbf8f93SLisandro Dalcin     if (stepok) break;
65cdbf8f93SLisandro Dalcin   }
66be5899b3SLisandro Dalcin   if (reject >= ts->max_reject) {
67be5899b3SLisandro Dalcin     ts->reason = TS_DIVERGED_STEP_REJECTED;
6863a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject));
693ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
707bf11e45SBarry Smith   }
71be5899b3SLisandro Dalcin 
72be5899b3SLisandro Dalcin   ts->ptime += ts->time_step;
73be5899b3SLisandro Dalcin   ts->time_step = next_time_step;
74be5899b3SLisandro Dalcin 
759566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(pseudo->xdot));
76eb636060SBarry Smith   PetscCall(TSComputeIFunction(ts, ts->ptime, pseudo->update, pseudo->xdot, pseudo->func, PETSC_FALSE));
77eb636060SBarry Smith   PetscCall(PetscObjectStateGet((PetscObject)pseudo->update, &pseudo->Xstate));
789566063dSJacob Faibussowitsch   PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
79eb636060SBarry Smith 
803118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
813118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
8263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol));
833ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
843118ae5eSBarry Smith   }
853118ae5eSBarry Smith   if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) {
863118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
8763a3b9bcSJacob 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));
883ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
893118ae5eSBarry Smith   }
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
912d3f70b5SBarry Smith }
922d3f70b5SBarry Smith 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts)
94d71ae5a4SJacob Faibussowitsch {
957bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
962d3f70b5SBarry Smith 
973a40ed3dSBarry Smith   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->update));
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->func));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->xdot));
1013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1022d3f70b5SBarry Smith }
1032d3f70b5SBarry Smith 
104d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts)
105d71ae5a4SJacob Faibussowitsch {
106277b19d0SLisandro Dalcin   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(TSReset_Pseudo(ts));
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL));
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL));
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115277b19d0SLisandro Dalcin }
1162d3f70b5SBarry Smith 
117d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot)
118d71ae5a4SJacob Faibussowitsch {
1196f2d6a7bSJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1202d3f70b5SBarry Smith 
1213a40ed3dSBarry Smith   PetscFunctionBegin;
1226f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1242d3f70b5SBarry Smith }
1252d3f70b5SBarry Smith 
1266f2d6a7bSJed Brown /*
1276f2d6a7bSJed Brown     The transient residual is
1286f2d6a7bSJed Brown 
1296f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
1306f2d6a7bSJed Brown 
1316f2d6a7bSJed Brown     or for ODE,
1326f2d6a7bSJed Brown 
1336f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
1346f2d6a7bSJed Brown 
1356f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
1366f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
1376f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
1386f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
1396f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
140eb636060SBarry Smith     algorithm, it only takes this one full Newton step with the steady state
1416f2d6a7bSJed Brown     residual, and then advances to the next time step.
142eb636060SBarry Smith 
143eb636060SBarry Smith     This routine avoids a repeated identical call to TSComputeRHSFunction() when that result
144eb636060SBarry Smith     is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo()
145eb636060SBarry Smith 
1466f2d6a7bSJed Brown */
147d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts)
148d71ae5a4SJacob Faibussowitsch {
149eb636060SBarry Smith   TSIFunctionFn    *ifunction;
150eb636060SBarry Smith   TSRHSFunctionFn  *rhsfunction;
151eb636060SBarry Smith   void             *ctx;
152eb636060SBarry Smith   DM                dm;
153eb636060SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo *)ts->data;
15481775c15SStefano Zampini   const PetscScalar mdt    = 1.0 / ts->time_step;
155eb636060SBarry Smith   PetscBool         KSPSNES;
156eb636060SBarry Smith   PetscObjectState  Xstate;
1572d3f70b5SBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
159eb636060SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
160eb636060SBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
161eb636060SBarry Smith   PetscValidHeaderSpecific(Y, VEC_CLASSID, 3);
162eb636060SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 4);
163eb636060SBarry Smith 
164eb636060SBarry Smith   PetscCall(TSGetDM(ts, &dm));
165eb636060SBarry Smith   PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx));
166eb636060SBarry Smith   PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL));
167eb636060SBarry Smith   PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()");
168eb636060SBarry Smith 
169eb636060SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESKSPONLY, &KSPSNES));
170eb636060SBarry Smith   PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate));
171eb636060SBarry Smith   if (Xstate != pseudo->Xstate || ifunction || !KSPSNES) {
172*c558785fSJames Wright     // X == ts->vec_sol for first SNES iteration, so pseudo->xdot == 0
17381775c15SStefano Zampini     PetscCall(VecAXPBYPCZ(pseudo->xdot, -mdt, mdt, 0, ts->vec_sol, X));
17481775c15SStefano Zampini     PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, pseudo->xdot, Y, PETSC_FALSE));
175eb636060SBarry Smith   } else {
176eb636060SBarry Smith     /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */
177eb636060SBarry Smith     /* note that pseudo->func contains the negation of TSComputeRHSFunction() */
17881775c15SStefano Zampini     PetscCall(VecWAXPY(Y, 1, pseudo->func, pseudo->xdot));
179eb636060SBarry Smith   }
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1816f2d6a7bSJed Brown }
1822d3f70b5SBarry Smith 
1836f2d6a7bSJed Brown /*
1846f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
1856f2d6a7bSJed Brown 
1866f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
1876f2d6a7bSJed Brown 
1886f2d6a7bSJed Brown     and for ODE:
1896f2d6a7bSJed Brown 
1906f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
1916f2d6a7bSJed Brown */
192d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts)
193d71ae5a4SJacob Faibussowitsch {
1946f2d6a7bSJed Brown   Vec Xdot;
1956f2d6a7bSJed Brown 
1966f2d6a7bSJed Brown   PetscFunctionBegin;
19781775c15SStefano Zampini   /* Xdot has already been computed in SNESTSFormFunction_Pseudo (SNES guarantees this) */
1989566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
1999566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2012d3f70b5SBarry Smith }
2022d3f70b5SBarry Smith 
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts)
204d71ae5a4SJacob Faibussowitsch {
2057bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
2062d3f70b5SBarry Smith 
2073a40ed3dSBarry Smith   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update));
2099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func));
2109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot));
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2122d3f70b5SBarry Smith }
21381775c15SStefano Zampini 
214d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy)
215d71ae5a4SJacob Faibussowitsch {
2167bf11e45SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
217ce94432eSBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
2182d3f70b5SBarry Smith 
2193a40ed3dSBarry Smith   PetscFunctionBegin;
220193ac0bcSJed Brown   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
2219566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
2229566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
223eb636060SBarry Smith     PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate));
2249566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
225193ac0bcSJed Brown   }
2269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
22763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm));
2289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2302d3f70b5SBarry Smith }
2312d3f70b5SBarry Smith 
232ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject)
233d71ae5a4SJacob Faibussowitsch {
2344bbc92c1SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
235ace3abfcSBarry Smith   PetscBool   flg    = PETSC_FALSE;
236649052a6SBarry Smith   PetscViewer viewer;
2372d3f70b5SBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
239d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options");
2409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL));
2412d3f70b5SBarry Smith   if (flg) {
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer));
24349abdd8aSBarry Smith     PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
24428aa8177SBarry Smith   }
245be5899b3SLisandro Dalcin   flg = pseudo->increment_dt_from_initial_dt;
2469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL));
247be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
2489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL));
2499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL));
2509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL));
2519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL));
252d0609cedSBarry Smith   PetscOptionsHeadEnd();
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2542d3f70b5SBarry Smith }
2552d3f70b5SBarry Smith 
256d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer)
257d71ae5a4SJacob Faibussowitsch {
2583118ae5eSBarry Smith   PetscBool isascii;
259d52bd9f3SBarry Smith 
2603a40ed3dSBarry Smith   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2623118ae5eSBarry Smith   if (isascii) {
2633118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
2649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  frtol - relative tolerance in function value %g\n", (double)pseudo->frtol));
2659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol));
2669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_initial - initial timestep %g\n", (double)pseudo->dt_initial));
2679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment));
2689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_max - maximum time %g\n", (double)pseudo->dt_max));
2693118ae5eSBarry Smith   }
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2712d3f70b5SBarry Smith }
2722d3f70b5SBarry Smith 
273ac226902SBarry Smith /*@C
274*c558785fSJames Wright   TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.  Use with `TSPseudoSetTimeStep()`.
275*c558785fSJames Wright 
276*c558785fSJames Wright   Collective, No Fortran Support
277*c558785fSJames Wright 
278*c558785fSJames Wright   Input Parameters:
279*c558785fSJames Wright + ts    - the timestep context
280*c558785fSJames Wright - dtctx - unused timestep context
281*c558785fSJames Wright 
282*c558785fSJames Wright   Output Parameter:
283*c558785fSJames Wright . newdt - the timestep to use for the next step
284*c558785fSJames Wright 
285*c558785fSJames Wright   Level: advanced
286*c558785fSJames Wright 
287*c558785fSJames Wright .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO`
288*c558785fSJames Wright @*/
289*c558785fSJames Wright PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx)
290*c558785fSJames Wright {
291*c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
292*c558785fSJames Wright   PetscReal  inc    = pseudo->dt_increment;
293*c558785fSJames Wright 
294*c558785fSJames Wright   PetscFunctionBegin;
295*c558785fSJames Wright   if (pseudo->fnorm < 0.0) {
296*c558785fSJames Wright     PetscCall(VecZeroEntries(pseudo->xdot));
297*c558785fSJames Wright     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
298*c558785fSJames Wright     PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate));
299*c558785fSJames Wright     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
300*c558785fSJames Wright   }
301*c558785fSJames Wright   if (pseudo->fnorm_initial < 0) {
302*c558785fSJames Wright     /* first time through so compute initial function norm */
303*c558785fSJames Wright     pseudo->fnorm_initial  = pseudo->fnorm;
304*c558785fSJames Wright     pseudo->fnorm_previous = pseudo->fnorm;
305*c558785fSJames Wright   }
306*c558785fSJames Wright   if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step;
307*c558785fSJames Wright   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm;
308*c558785fSJames Wright   else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm;
309*c558785fSJames Wright   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max);
310*c558785fSJames Wright   pseudo->fnorm_previous = pseudo->fnorm;
311*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
312*c558785fSJames Wright }
313*c558785fSJames Wright 
314*c558785fSJames 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)
315*c558785fSJames Wright {
316*c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
317*c558785fSJames Wright 
318*c558785fSJames Wright   PetscFunctionBegin;
319*c558785fSJames Wright   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
320*c558785fSJames Wright   PetscCall((*pseudo->dt)(ts, next_h, pseudo->dtctx));
321*c558785fSJames Wright   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
322*c558785fSJames Wright 
323*c558785fSJames Wright   *next_sc = 0;
324*c558785fSJames Wright   *wlte    = -1; /* Weighted local truncation error was not evaluated */
325*c558785fSJames Wright   *wltea   = -1; /* Weighted absolute local truncation error was not evaluated */
326*c558785fSJames Wright   *wlter   = -1; /* Weighted relative local truncation error was not evaluated */
327*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
328*c558785fSJames Wright }
329*c558785fSJames Wright 
330*c558785fSJames Wright static PetscErrorCode TSAdaptCheckStage_TSPseudo(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
331*c558785fSJames Wright {
332*c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
333*c558785fSJames Wright 
334*c558785fSJames Wright   PetscFunctionBegin;
335*c558785fSJames Wright   if (pseudo->verify) {
336*c558785fSJames Wright     PetscReal dt;
337*c558785fSJames Wright     PetscCall(TSGetTimeStep(ts, &dt));
338*c558785fSJames Wright     PetscCall((*pseudo->verify)(ts, Y, pseudo->verifyctx, &dt, accept));
339*c558785fSJames Wright     PetscCall(TSSetTimeStep(ts, dt));
340*c558785fSJames Wright   }
341*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
342*c558785fSJames Wright }
343*c558785fSJames Wright 
344*c558785fSJames Wright /*MC
345*c558785fSJames Wright   TSADAPTTSPSEUDO - TSPseudo adaptive controller for time stepping
346*c558785fSJames Wright 
347*c558785fSJames Wright   Level: developer
348*c558785fSJames Wright 
349*c558785fSJames Wright   Note:
350*c558785fSJames Wright   This is only meant to be used with `TSPSEUDO` time integrator.
351*c558785fSJames Wright 
352*c558785fSJames Wright .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`, `TSPSEUDO`
353*c558785fSJames Wright M*/
354*c558785fSJames Wright static PetscErrorCode TSAdaptCreate_TSPseudo(TSAdapt adapt)
355*c558785fSJames Wright {
356*c558785fSJames Wright   PetscFunctionBegin;
357*c558785fSJames Wright   adapt->ops->choose = TSAdaptChoose_TSPseudo;
358*c558785fSJames Wright   adapt->checkstage  = TSAdaptCheckStage_TSPseudo;
359*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
360*c558785fSJames Wright }
361*c558785fSJames Wright 
362*c558785fSJames Wright /*@
363*c558785fSJames Wright   TSPseudoComputeTimeStep - Computes the next timestep for a currently running
364*c558785fSJames Wright   pseudo-timestepping process.
365*c558785fSJames Wright 
366*c558785fSJames Wright   Collective
367*c558785fSJames Wright 
368*c558785fSJames Wright   Input Parameter:
369*c558785fSJames Wright . ts - timestep context
370*c558785fSJames Wright 
371*c558785fSJames Wright   Output Parameter:
372*c558785fSJames Wright . dt - newly computed timestep
373*c558785fSJames Wright 
374*c558785fSJames Wright   Level: developer
375*c558785fSJames Wright 
376*c558785fSJames Wright   Note:
377*c558785fSJames Wright   The routine to be called here to compute the timestep should be
378*c558785fSJames Wright   set by calling `TSPseudoSetTimeStep()`.
379*c558785fSJames Wright 
380*c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
381*c558785fSJames Wright @*/
382*c558785fSJames Wright PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt)
383*c558785fSJames Wright {
384*c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
385*c558785fSJames Wright 
386*c558785fSJames Wright   PetscFunctionBegin;
387*c558785fSJames Wright   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
388*c558785fSJames Wright   PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx));
389*c558785fSJames Wright   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
390*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
391*c558785fSJames Wright }
392*c558785fSJames Wright 
393*c558785fSJames Wright /*@C
394*c558785fSJames Wright   TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
395*c558785fSJames Wright 
396*c558785fSJames Wright   Collective, No Fortran Support
397*c558785fSJames Wright 
398*c558785fSJames Wright   Input Parameters:
399*c558785fSJames Wright + ts     - the timestep context
400*c558785fSJames Wright . dtctx  - unused timestep context
401*c558785fSJames Wright - update - latest solution vector
402*c558785fSJames Wright 
403*c558785fSJames Wright   Output Parameters:
404*c558785fSJames Wright + newdt - the timestep to use for the next step
405*c558785fSJames Wright - flag  - flag indicating whether the last time step was acceptable
406*c558785fSJames Wright 
407*c558785fSJames Wright   Level: advanced
408*c558785fSJames Wright 
409*c558785fSJames Wright   Note:
410*c558785fSJames Wright   This routine always returns a flag of 1, indicating an acceptable
411*c558785fSJames Wright   timestep.
412*c558785fSJames Wright 
413*c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
414*c558785fSJames Wright @*/
415*c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag)
416*c558785fSJames Wright {
417*c558785fSJames Wright   PetscFunctionBegin;
418*c558785fSJames Wright   // NOTE: This function was never used
419*c558785fSJames Wright   *flag = PETSC_TRUE;
420*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
421*c558785fSJames Wright }
422*c558785fSJames Wright 
423*c558785fSJames Wright /*@
424*c558785fSJames Wright   TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
425*c558785fSJames Wright 
426*c558785fSJames Wright   Collective
427*c558785fSJames Wright 
428*c558785fSJames Wright   Input Parameters:
429*c558785fSJames Wright + ts     - timestep context
430*c558785fSJames Wright - update - latest solution vector
431*c558785fSJames Wright 
432*c558785fSJames Wright   Output Parameters:
433*c558785fSJames Wright + dt   - newly computed timestep (if it had to shrink)
434*c558785fSJames Wright - flag - indicates if current timestep was ok
435*c558785fSJames Wright 
436*c558785fSJames Wright   Level: advanced
437*c558785fSJames Wright 
438*c558785fSJames Wright   Notes:
439*c558785fSJames Wright   The routine to be called here to compute the timestep should be
440*c558785fSJames Wright   set by calling `TSPseudoSetVerifyTimeStep()`.
441*c558785fSJames Wright 
442*c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
443*c558785fSJames Wright @*/
444*c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag)
445*c558785fSJames Wright {
446*c558785fSJames Wright   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
447*c558785fSJames Wright 
448*c558785fSJames Wright   PetscFunctionBegin;
449*c558785fSJames Wright   // NOTE: This function is never used
450*c558785fSJames Wright   *flag = PETSC_TRUE;
451*c558785fSJames Wright   if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag));
452*c558785fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
453*c558785fSJames Wright }
454*c558785fSJames Wright 
455*c558785fSJames Wright /*@C
45682bf6240SBarry Smith   TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
45782bf6240SBarry Smith   last timestep.
45882bf6240SBarry Smith 
459c3339decSBarry Smith   Logically Collective
46015091d37SBarry Smith 
46182bf6240SBarry Smith   Input Parameters:
46215091d37SBarry Smith + ts  - timestep context
46382bf6240SBarry Smith . dt  - user-defined function to verify timestep
4642fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`)
46582bf6240SBarry Smith 
46620f4b53cSBarry Smith   Calling sequence of `func`:
4672fe279fdSBarry Smith + ts     - the time-step context
4682fe279fdSBarry Smith . update - latest solution vector
46914d0ab18SJacob Faibussowitsch . ctx    - [optional] user-defined timestep context
47082bf6240SBarry Smith . newdt  - the timestep to use for the next step
471a2b725a8SWilliam Gropp - flag   - flag indicating whether the last time step was acceptable
47282bf6240SBarry Smith 
473bcf0153eSBarry Smith   Level: advanced
474bcf0153eSBarry Smith 
475bcf0153eSBarry Smith   Note:
476bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoVerifyTimeStep()`
47782bf6240SBarry Smith   during the timestepping process.
47882bf6240SBarry Smith 
4791cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
48082bf6240SBarry Smith @*/
48114d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx)
482d71ae5a4SJacob Faibussowitsch {
48382bf6240SBarry Smith   PetscFunctionBegin;
4840700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
485cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx));
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48782bf6240SBarry Smith }
48882bf6240SBarry Smith 
48982bf6240SBarry Smith /*@
49082bf6240SBarry Smith   TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
4918d359177SBarry Smith   dt when using the TSPseudoTimeStepDefault() routine.
49282bf6240SBarry Smith 
493c3339decSBarry Smith   Logically Collective
494fee21e36SBarry Smith 
49515091d37SBarry Smith   Input Parameters:
49615091d37SBarry Smith + ts  - the timestep context
49715091d37SBarry Smith - inc - the scaling factor >= 1.0
49815091d37SBarry Smith 
49982bf6240SBarry Smith   Options Database Key:
50067b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment
50182bf6240SBarry Smith 
50215091d37SBarry Smith   Level: advanced
50315091d37SBarry Smith 
5041cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
50582bf6240SBarry Smith @*/
506d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc)
507d71ae5a4SJacob Faibussowitsch {
50882bf6240SBarry Smith   PetscFunctionBegin;
5090700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
510c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts, inc, 2);
511cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc));
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51382bf6240SBarry Smith }
51482bf6240SBarry Smith 
51586534af1SJed Brown /*@
51686534af1SJed Brown   TSPseudoSetMaxTimeStep - Sets the maximum time step
5178d359177SBarry Smith   when using the TSPseudoTimeStepDefault() routine.
51886534af1SJed Brown 
519c3339decSBarry Smith   Logically Collective
52086534af1SJed Brown 
52186534af1SJed Brown   Input Parameters:
52286534af1SJed Brown + ts    - the timestep context
52386534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate
52486534af1SJed Brown 
52586534af1SJed Brown   Options Database Key:
52667b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt
52786534af1SJed Brown 
52886534af1SJed Brown   Level: advanced
52986534af1SJed Brown 
5301cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
53186534af1SJed Brown @*/
532d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt)
533d71ae5a4SJacob Faibussowitsch {
53486534af1SJed Brown   PetscFunctionBegin;
53586534af1SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
53686534af1SJed Brown   PetscValidLogicalCollectiveReal(ts, maxdt, 2);
537cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
5383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53986534af1SJed Brown }
54086534af1SJed Brown 
54182bf6240SBarry Smith /*@
54282bf6240SBarry Smith   TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
543b44f4de4SBarry 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.$
54482bf6240SBarry Smith 
545c3339decSBarry Smith   Logically Collective
54615091d37SBarry Smith 
54782bf6240SBarry Smith   Input Parameter:
54882bf6240SBarry Smith . ts - the timestep context
54982bf6240SBarry Smith 
55082bf6240SBarry Smith   Options Database Key:
551b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment
55282bf6240SBarry Smith 
55315091d37SBarry Smith   Level: advanced
55415091d37SBarry Smith 
5551cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
55682bf6240SBarry Smith @*/
557d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
558d71ae5a4SJacob Faibussowitsch {
55982bf6240SBarry Smith   PetscFunctionBegin;
5600700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
561cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts));
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56382bf6240SBarry Smith }
56482bf6240SBarry Smith 
565ac226902SBarry Smith /*@C
56682bf6240SBarry Smith   TSPseudoSetTimeStep - Sets the user-defined routine to be
56782bf6240SBarry Smith   called at each pseudo-timestep to update the timestep.
56882bf6240SBarry Smith 
569c3339decSBarry Smith   Logically Collective
57015091d37SBarry Smith 
57182bf6240SBarry Smith   Input Parameters:
57215091d37SBarry Smith + ts  - timestep context
57382bf6240SBarry Smith . dt  - function to compute timestep
5742fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`)
57582bf6240SBarry Smith 
57620f4b53cSBarry Smith   Calling sequence of `dt`:
57714d0ab18SJacob Faibussowitsch + ts    - the `TS` context
57814d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep
57914d0ab18SJacob Faibussowitsch - ctx   - [optional] user-defined context
58082bf6240SBarry Smith 
581bcf0153eSBarry Smith   Level: intermediate
58282bf6240SBarry Smith 
583bcf0153eSBarry Smith   Notes:
584bcf0153eSBarry Smith   The routine set here will be called by `TSPseudoComputeTimeStep()`
585bcf0153eSBarry Smith   during the timestepping process.
586bcf0153eSBarry Smith 
587bcf0153eSBarry Smith   If not set then `TSPseudoTimeStepDefault()` is automatically used
588bcf0153eSBarry Smith 
5891cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
59082bf6240SBarry Smith @*/
59114d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx)
592d71ae5a4SJacob Faibussowitsch {
59382bf6240SBarry Smith   PetscFunctionBegin;
5940700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
595cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx));
5963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59782bf6240SBarry Smith }
59882bf6240SBarry Smith 
599ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/
600d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx)
601d71ae5a4SJacob Faibussowitsch {
602be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
60382bf6240SBarry Smith 
60482bf6240SBarry Smith   PetscFunctionBegin;
60582bf6240SBarry Smith   pseudo->verify    = dt;
60682bf6240SBarry Smith   pseudo->verifyctx = ctx;
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60882bf6240SBarry Smith }
60982bf6240SBarry Smith 
610d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc)
611d71ae5a4SJacob Faibussowitsch {
6124bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
61382bf6240SBarry Smith 
61482bf6240SBarry Smith   PetscFunctionBegin;
61582bf6240SBarry Smith   pseudo->dt_increment = inc;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61782bf6240SBarry Smith }
61882bf6240SBarry Smith 
619d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt)
620d71ae5a4SJacob Faibussowitsch {
62186534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
62286534af1SJed Brown 
62386534af1SJed Brown   PetscFunctionBegin;
62486534af1SJed Brown   pseudo->dt_max = maxdt;
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62686534af1SJed Brown }
62786534af1SJed Brown 
628d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
629d71ae5a4SJacob Faibussowitsch {
6304bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
63182bf6240SBarry Smith 
63282bf6240SBarry Smith   PetscFunctionBegin;
6334bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
6343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63582bf6240SBarry Smith }
63682bf6240SBarry Smith 
6376849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/
638d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx)
639d71ae5a4SJacob Faibussowitsch {
6404bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
64182bf6240SBarry Smith 
64282bf6240SBarry Smith   PetscFunctionBegin;
64382bf6240SBarry Smith   pseudo->dt    = dt;
64482bf6240SBarry Smith   pseudo->dtctx = ctx;
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64682bf6240SBarry Smith }
64782bf6240SBarry Smith 
64810e6a065SJed Brown /*MC
6491d27aa22SBarry Smith   TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97`
65082bf6240SBarry Smith 
65110e6a065SJed Brown   This method solves equations of the form
65210e6a065SJed Brown 
6531d27aa22SBarry Smith   $$
6541d27aa22SBarry Smith   F(X,Xdot) = 0
6551d27aa22SBarry Smith   $$
65610e6a065SJed Brown 
65710e6a065SJed Brown   for steady state using the iteration
65810e6a065SJed Brown 
6591d27aa22SBarry Smith   $$
66037fdd005SBarry Smith   [G'] S = -F(X,0)
66137fdd005SBarry Smith   X += S
6621d27aa22SBarry Smith   $$
66310e6a065SJed Brown 
66410e6a065SJed Brown   where
66510e6a065SJed Brown 
6661d27aa22SBarry Smith   $$
6671d27aa22SBarry Smith   G(Y) = F(Y,(Y-X)/dt)
6681d27aa22SBarry Smith   $$
66910e6a065SJed Brown 
6706f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6716f2d6a7bSJed Brown   state".  See note below.
67210e6a065SJed Brown 
67393a54799SPierre Jolivet   In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`.
674fb59f417SJames Wright   It determines the next timestep via
675fb59f417SJames Wright 
676fb59f417SJames Wright   $$
677fb59f417SJames Wright   dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert}
678fb59f417SJames Wright   $$
679fb59f417SJames Wright 
680fb59f417SJames Wright   where $r$ is an additional growth factor (set by `-ts_pseudo_increment`).
681fb59f417SJames Wright   An alternative formulation is also available that uses the initial timestep and function norm.
682fb59f417SJames Wright 
683fb59f417SJames Wright   $$
684fb59f417SJames Wright   dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert}
685fb59f417SJames Wright   $$
686fb59f417SJames Wright 
687fb59f417SJames Wright   This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`.
688fb59f417SJames Wright   For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`.
689fb59f417SJames Wright 
690bcf0153eSBarry Smith   Options Database Keys:
69110e6a065SJed Brown +  -ts_pseudo_increment <real>                     - ratio of increase dt
6923118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
693fb59f417SJames Wright .  -ts_pseudo_max_dt                               - Maximum dt for adaptive timestepping algorithm
694fb59f417SJames Wright .  -ts_pseudo_monitor                              - Monitor convergence of the function norm
695fb59f417SJames Wright .  -ts_pseudo_fatol <atol>                         - stop iterating when the function norm is less than `atol`
696fb59f417SJames Wright -  -ts_pseudo_frtol <rtol>                         - stop iterating when the function norm divided by the initial function norm is less than `rtol`
69710e6a065SJed Brown 
69810e6a065SJed Brown   Level: beginner
69910e6a065SJed Brown 
70010e6a065SJed Brown   Notes:
7016f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
7026f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
7036f2d6a7bSJed Brown   last step, on the first Newton iteration we have
7046f2d6a7bSJed Brown 
7051d27aa22SBarry Smith   $$
7061d27aa22SBarry Smith   Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
7071d27aa22SBarry Smith   $$
7086f2d6a7bSJed Brown 
709eb636060SBarry Smith   The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it
710eb636060SBarry Smith   is $ dF/dX + shift*1/dt $. Hence  still contains the $ 1/dt $ term so the Jacobian is not simply the
711eb636060SBarry Smith   Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$.
712eb636060SBarry Smith 
7136f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
7146f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
7156f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
716fb59f417SJames Wright   By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers.
717fb59f417SJames Wright   Setting the `SNESType` via `-snes_type` will override this default setting.
71810e6a065SJed Brown 
7191cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`
72010e6a065SJed Brown M*/
721d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
722d71ae5a4SJacob Faibussowitsch {
7237bf11e45SBarry Smith   TS_Pseudo *pseudo;
724193ac0bcSJed Brown   SNES       snes;
72519fd82e9SBarry Smith   SNESType   stype;
7262d3f70b5SBarry Smith 
7273a40ed3dSBarry Smith   PetscFunctionBegin;
728277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
729000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
730000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
731000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
732000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
733000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
7340f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
7350f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
736*c558785fSJames Wright   ts->default_adapt_type  = TSADAPTTSPSEUDO;
737825ab935SBarry Smith   ts->usessnes            = PETSC_TRUE;
7387bf11e45SBarry Smith 
739*c558785fSJames Wright   PetscCall(TSAdaptRegister(TSADAPTTSPSEUDO, TSAdaptCreate_TSPseudo));
740*c558785fSJames Wright 
7419566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
7429566063dSJacob Faibussowitsch   PetscCall(SNESGetType(snes, &stype));
7439566063dSJacob Faibussowitsch   if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY));
7442d3f70b5SBarry Smith 
7454dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pseudo));
7467bf11e45SBarry Smith   ts->data = (void *)pseudo;
7472d3f70b5SBarry Smith 
748be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
749be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
75028aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7514bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
752193ac0bcSJed Brown   pseudo->fnorm                        = -1;
753be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
754be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
7553118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7563118ae5eSBarry Smith   pseudo->fatol = 1.e-25;
7573118ae5eSBarry Smith   pseudo->frtol = 1.e-5;
7583118ae5eSBarry Smith #else
7593118ae5eSBarry Smith   pseudo->fatol = 1.e-50;
7603118ae5eSBarry Smith   pseudo->frtol = 1.e-12;
7613118ae5eSBarry Smith #endif
7629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo));
7639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo));
7649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo));
7659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo));
7669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo));
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7682d3f70b5SBarry Smith }
769