xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
52d3f70b5SBarry Smith 
62d3f70b5SBarry Smith typedef struct {
72d3f70b5SBarry Smith   Vec update; /* work vector where new solution is formed */
82d3f70b5SBarry Smith   Vec func;   /* work vector where F(t[i],u[i]) is stored */
96f2d6a7bSJed Brown   Vec xdot;   /* work vector for time derivative of state */
102d3f70b5SBarry Smith 
112d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
122d3f70b5SBarry Smith 
136849ba73SBarry Smith   PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */
142d3f70b5SBarry Smith   void *dtctx;
15ace3abfcSBarry Smith   PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */
167bf11e45SBarry Smith   void *verifyctx;
172d3f70b5SBarry Smith 
18cdbf8f93SLisandro Dalcin   PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */
1987828ca2SBarry Smith   PetscReal fnorm_previous;
2028aa8177SBarry Smith 
21cdbf8f93SLisandro Dalcin   PetscReal dt_initial;   /* initial time-step */
2287828ca2SBarry Smith   PetscReal dt_increment; /* scaling that dt is incremented each time-step */
2386534af1SJed Brown   PetscReal dt_max;       /* maximum time step */
24ace3abfcSBarry Smith   PetscBool increment_dt_from_initial_dt;
253118ae5eSBarry Smith   PetscReal fatol, frtol;
267bf11e45SBarry Smith } TS_Pseudo;
272d3f70b5SBarry Smith 
282d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
292d3f70b5SBarry Smith 
308d359177SBarry Smith /*@C
317bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
332d3f70b5SBarry Smith 
3415091d37SBarry Smith     Collective on TS
3515091d37SBarry Smith 
367bf11e45SBarry Smith     Input Parameter:
377bf11e45SBarry Smith .   ts - timestep context
387bf11e45SBarry Smith 
397bf11e45SBarry Smith     Output Parameter:
40fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
41fb4a63b6SLois Curfman McInnes 
428d359177SBarry Smith     Level: developer
43564e8f4eSLois Curfman McInnes 
44564e8f4eSLois Curfman McInnes     Notes:
45564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
46564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
47564e8f4eSLois Curfman McInnes 
48db781477SPatrick Sanan .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
497bf11e45SBarry Smith @*/
509371c9d4SSatish Balay PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) {
517bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
527bf11e45SBarry Smith 
533a40ed3dSBarry Smith   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
559566063dSJacob Faibussowitsch   PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx));
569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0));
573a40ed3dSBarry Smith   PetscFunctionReturn(0);
587bf11e45SBarry Smith }
597bf11e45SBarry Smith 
607bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
617bf11e45SBarry Smith /*@C
628d359177SBarry Smith    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
637bf11e45SBarry Smith 
6415091d37SBarry Smith    Collective on TS
6515091d37SBarry Smith 
667bf11e45SBarry Smith    Input Parameters:
6715091d37SBarry Smith +  ts - the timestep context
687bf11e45SBarry Smith .  dtctx - unused timestep context
6915091d37SBarry Smith -  update - latest solution vector
707bf11e45SBarry Smith 
71564e8f4eSLois Curfman McInnes    Output Parameters:
7215091d37SBarry Smith +  newdt - the timestep to use for the next step
7315091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
747bf11e45SBarry Smith 
7515091d37SBarry Smith    Level: advanced
76fee21e36SBarry Smith 
77564e8f4eSLois Curfman McInnes    Note:
78564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
79564e8f4eSLois Curfman McInnes    timestep.
80564e8f4eSLois Curfman McInnes 
81db781477SPatrick Sanan .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
827bf11e45SBarry Smith @*/
839371c9d4SSatish Balay PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) {
843a40ed3dSBarry Smith   PetscFunctionBegin;
85a7cc72afSBarry Smith   *flag = PETSC_TRUE;
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
877bf11e45SBarry Smith }
887bf11e45SBarry Smith 
897bf11e45SBarry Smith /*@
90564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
917bf11e45SBarry Smith 
9215091d37SBarry Smith     Collective on TS
9315091d37SBarry Smith 
94fb4a63b6SLois Curfman McInnes     Input Parameters:
9515091d37SBarry Smith +   ts - timestep context
9615091d37SBarry Smith -   update - latest solution vector
977bf11e45SBarry Smith 
98fb4a63b6SLois Curfman McInnes     Output Parameters:
9915091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
10015091d37SBarry Smith -   flag - indicates if current timestep was ok
1017bf11e45SBarry Smith 
10215091d37SBarry Smith     Level: advanced
103fee21e36SBarry Smith 
104564e8f4eSLois Curfman McInnes     Notes:
105564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
106564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
107564e8f4eSLois Curfman McInnes 
108db781477SPatrick Sanan .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
1097bf11e45SBarry Smith @*/
1109371c9d4SSatish Balay PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) {
1117bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1127bf11e45SBarry Smith 
1133a40ed3dSBarry Smith   PetscFunctionBegin;
114cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
1151baa6e33SBarry Smith   if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag));
1163a40ed3dSBarry Smith   PetscFunctionReturn(0);
1177bf11e45SBarry Smith }
1187bf11e45SBarry Smith 
1197bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1207bf11e45SBarry Smith 
1219371c9d4SSatish Balay static PetscErrorCode TSStep_Pseudo(TS ts) {
122277b19d0SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
123be5899b3SLisandro Dalcin   PetscInt   nits, lits, reject;
124cdbf8f93SLisandro Dalcin   PetscBool  stepok;
125be5899b3SLisandro Dalcin   PetscReal  next_time_step = ts->time_step;
1262d3f70b5SBarry Smith 
1273a40ed3dSBarry Smith   PetscFunctionBegin;
128bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
1299566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, pseudo->update));
1309566063dSJacob Faibussowitsch   PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step));
131cdbf8f93SLisandro Dalcin   for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) {
132cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
1339566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime + ts->time_step));
1349566063dSJacob Faibussowitsch     PetscCall(SNESSolve(ts->snes, NULL, pseudo->update));
1359566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1369566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1379371c9d4SSatish Balay     ts->snes_its += nits;
1389371c9d4SSatish Balay     ts->ksp_its += lits;
1399566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &(pseudo->update)));
1409566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok));
1419371c9d4SSatish Balay     if (!stepok) {
1429371c9d4SSatish Balay       next_time_step = ts->time_step;
1439371c9d4SSatish Balay       continue;
1449371c9d4SSatish Balay     }
145193ac0bcSJed Brown     pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
1469566063dSJacob Faibussowitsch     PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok));
147cdbf8f93SLisandro Dalcin     if (stepok) break;
148cdbf8f93SLisandro Dalcin   }
149be5899b3SLisandro Dalcin   if (reject >= ts->max_reject) {
150be5899b3SLisandro Dalcin     ts->reason = TS_DIVERGED_STEP_REJECTED;
15163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject));
152cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1537bf11e45SBarry Smith   }
154be5899b3SLisandro Dalcin 
1559566063dSJacob Faibussowitsch   PetscCall(VecCopy(pseudo->update, ts->vec_sol));
156be5899b3SLisandro Dalcin   ts->ptime += ts->time_step;
157be5899b3SLisandro Dalcin   ts->time_step = next_time_step;
158be5899b3SLisandro Dalcin 
1593118ae5eSBarry Smith   if (pseudo->fnorm < 0) {
1609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
1619566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
1629566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
1633118ae5eSBarry Smith   }
1643118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1653118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
16663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol));
1673118ae5eSBarry Smith     PetscFunctionReturn(0);
1683118ae5eSBarry Smith   }
1693118ae5eSBarry Smith   if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) {
1703118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
17163a3b9bcSJacob 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));
1723118ae5eSBarry Smith     PetscFunctionReturn(0);
1733118ae5eSBarry Smith   }
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
1752d3f70b5SBarry Smith }
1762d3f70b5SBarry Smith 
1772d3f70b5SBarry Smith /*------------------------------------------------------------*/
1789371c9d4SSatish Balay static PetscErrorCode TSReset_Pseudo(TS ts) {
1797bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
1802d3f70b5SBarry Smith 
1813a40ed3dSBarry Smith   PetscFunctionBegin;
1829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->update));
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->func));
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->xdot));
1853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1862d3f70b5SBarry Smith }
1872d3f70b5SBarry Smith 
1889371c9d4SSatish Balay static PetscErrorCode TSDestroy_Pseudo(TS ts) {
189277b19d0SLisandro Dalcin   PetscFunctionBegin;
1909566063dSJacob Faibussowitsch   PetscCall(TSReset_Pseudo(ts));
1919566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
1929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL));
1939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL));
1949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL));
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL));
1969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL));
197277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
198277b19d0SLisandro Dalcin }
1992d3f70b5SBarry Smith 
2002d3f70b5SBarry Smith /*------------------------------------------------------------*/
2012d3f70b5SBarry Smith 
2026f2d6a7bSJed Brown /*
2036f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2046f2d6a7bSJed Brown */
2059371c9d4SSatish Balay static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) {
2066f2d6a7bSJed Brown   TS_Pseudo        *pseudo = (TS_Pseudo *)ts->data;
207193ac0bcSJed Brown   const PetscScalar mdt    = 1.0 / ts->time_step, *xnp1, *xn;
208193ac0bcSJed Brown   PetscScalar      *xdot;
209a7cc72afSBarry Smith   PetscInt          i, n;
2102d3f70b5SBarry Smith 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
212aab5bcd8SJed Brown   *Xdot = NULL;
2139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ts->vec_sol, &xn));
2149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xnp1));
2159566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pseudo->xdot, &xdot));
2169566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
217bbd56ea5SKarl Rupp   for (i = 0; i < n; i++) xdot[i] = mdt * (xnp1[i] - xn[i]);
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ts->vec_sol, &xn));
2199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xnp1));
2209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pseudo->xdot, &xdot));
2216f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2223a40ed3dSBarry Smith   PetscFunctionReturn(0);
2232d3f70b5SBarry Smith }
2242d3f70b5SBarry Smith 
2256f2d6a7bSJed Brown /*
2266f2d6a7bSJed Brown     The transient residual is
2276f2d6a7bSJed Brown 
2286f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2296f2d6a7bSJed Brown 
2306f2d6a7bSJed Brown     or for ODE,
2316f2d6a7bSJed Brown 
2326f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2336f2d6a7bSJed Brown 
2346f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2356f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2366f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2376f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2386f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2396f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2406f2d6a7bSJed Brown     residual, and then advances to the next time step.
2416f2d6a7bSJed Brown */
2429371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) {
2436f2d6a7bSJed Brown   Vec Xdot;
2442d3f70b5SBarry Smith 
2453a40ed3dSBarry Smith   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
2479566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, Xdot, Y, PETSC_FALSE));
2486f2d6a7bSJed Brown   PetscFunctionReturn(0);
2496f2d6a7bSJed Brown }
2502d3f70b5SBarry Smith 
2516f2d6a7bSJed Brown /*
2526f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2536f2d6a7bSJed Brown 
2546f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2556f2d6a7bSJed Brown 
2566f2d6a7bSJed Brown     and for ODE:
2576f2d6a7bSJed Brown 
2586f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2596f2d6a7bSJed Brown */
2609371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) {
2616f2d6a7bSJed Brown   Vec Xdot;
2626f2d6a7bSJed Brown 
2636f2d6a7bSJed Brown   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts, X, &Xdot));
2659566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE));
2663a40ed3dSBarry Smith   PetscFunctionReturn(0);
2672d3f70b5SBarry Smith }
2682d3f70b5SBarry Smith 
2699371c9d4SSatish Balay static PetscErrorCode TSSetUp_Pseudo(TS ts) {
2707bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
2712d3f70b5SBarry Smith 
2723a40ed3dSBarry Smith   PetscFunctionBegin;
2739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update));
2749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func));
2759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot));
2763a40ed3dSBarry Smith   PetscFunctionReturn(0);
2772d3f70b5SBarry Smith }
2782d3f70b5SBarry Smith /*------------------------------------------------------------*/
2792d3f70b5SBarry Smith 
2809371c9d4SSatish Balay static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) {
2817bf11e45SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
282ce94432eSBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
2832d3f70b5SBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
285193ac0bcSJed Brown   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
2869566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
2879566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
2889566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
289193ac0bcSJed Brown   }
2909566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
29163a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm));
2929566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
2933a40ed3dSBarry Smith   PetscFunctionReturn(0);
2942d3f70b5SBarry Smith }
2952d3f70b5SBarry Smith 
2969371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems *PetscOptionsObject) {
2974bbc92c1SBarry Smith   TS_Pseudo  *pseudo = (TS_Pseudo *)ts->data;
298ace3abfcSBarry Smith   PetscBool   flg    = PETSC_FALSE;
299649052a6SBarry Smith   PetscViewer viewer;
3002d3f70b5SBarry Smith 
3013a40ed3dSBarry Smith   PetscFunctionBegin;
302d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options");
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL));
3042d3f70b5SBarry Smith   if (flg) {
3059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer));
3069566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
30728aa8177SBarry Smith   }
308be5899b3SLisandro Dalcin   flg = pseudo->increment_dt_from_initial_dt;
3099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL));
310be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
3119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL));
3129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL));
3139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL));
3149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL));
315d0609cedSBarry Smith   PetscOptionsHeadEnd();
3163a40ed3dSBarry Smith   PetscFunctionReturn(0);
3172d3f70b5SBarry Smith }
3182d3f70b5SBarry Smith 
3199371c9d4SSatish Balay static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) {
3203118ae5eSBarry Smith   PetscBool isascii;
321d52bd9f3SBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
3239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3243118ae5eSBarry Smith   if (isascii) {
3253118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
3269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  frtol - relative tolerance in function value %g\n", (double)pseudo->frtol));
3279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol));
3289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_initial - initial timestep %g\n", (double)pseudo->dt_initial));
3299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment));
3309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dt_max - maximum time %g\n", (double)pseudo->dt_max));
3313118ae5eSBarry Smith   }
3323a40ed3dSBarry Smith   PetscFunctionReturn(0);
3332d3f70b5SBarry Smith }
3342d3f70b5SBarry Smith 
33582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
336ac226902SBarry Smith /*@C
33782bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
33882bf6240SBarry Smith    last timestep.
33982bf6240SBarry Smith 
3403f9fe445SBarry Smith    Logically Collective on TS
34115091d37SBarry Smith 
34282bf6240SBarry Smith    Input Parameters:
34315091d37SBarry Smith +  ts - timestep context
34482bf6240SBarry Smith .  dt - user-defined function to verify timestep
34515091d37SBarry Smith -  ctx - [optional] user-defined context for private data
3460298fd71SBarry Smith          for the timestep verification routine (may be NULL)
34782bf6240SBarry Smith 
34815091d37SBarry Smith    Level: advanced
349fee21e36SBarry Smith 
35082bf6240SBarry Smith    Calling sequence of func:
351a2b725a8SWilliam Gropp $  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
35282bf6240SBarry Smith 
353a2b725a8SWilliam Gropp +  update - latest solution vector
35482bf6240SBarry Smith .  ctx - [optional] timestep context
35582bf6240SBarry Smith .  newdt - the timestep to use for the next step
356a2b725a8SWilliam Gropp -  flag - flag indicating whether the last time step was acceptable
35782bf6240SBarry Smith 
35882bf6240SBarry Smith    Notes:
35982bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
36082bf6240SBarry Smith    during the timestepping process.
36182bf6240SBarry Smith 
362db781477SPatrick Sanan .seealso: `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
36382bf6240SBarry Smith @*/
3649371c9d4SSatish Balay PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS, Vec, void *, PetscReal *, PetscBool *), void *ctx) {
36582bf6240SBarry Smith   PetscFunctionBegin;
3660700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
367cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode(*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx));
36882bf6240SBarry Smith   PetscFunctionReturn(0);
36982bf6240SBarry Smith }
37082bf6240SBarry Smith 
37182bf6240SBarry Smith /*@
37282bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
3738d359177SBarry Smith     dt when using the TSPseudoTimeStepDefault() routine.
37482bf6240SBarry Smith 
3753f9fe445SBarry Smith    Logically Collective on TS
376fee21e36SBarry Smith 
37715091d37SBarry Smith     Input Parameters:
37815091d37SBarry Smith +   ts - the timestep context
37915091d37SBarry Smith -   inc - the scaling factor >= 1.0
38015091d37SBarry Smith 
38182bf6240SBarry Smith     Options Database Key:
38267b8a455SSatish Balay .    -ts_pseudo_increment <increment> - set pseudo increment
38382bf6240SBarry Smith 
38415091d37SBarry Smith     Level: advanced
38515091d37SBarry Smith 
386db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
38782bf6240SBarry Smith @*/
3889371c9d4SSatish Balay PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) {
38982bf6240SBarry Smith   PetscFunctionBegin;
3900700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
391c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts, inc, 2);
392cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc));
39382bf6240SBarry Smith   PetscFunctionReturn(0);
39482bf6240SBarry Smith }
39582bf6240SBarry Smith 
39686534af1SJed Brown /*@
39786534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
3988d359177SBarry Smith     when using the TSPseudoTimeStepDefault() routine.
39986534af1SJed Brown 
40086534af1SJed Brown    Logically Collective on TS
40186534af1SJed Brown 
40286534af1SJed Brown     Input Parameters:
40386534af1SJed Brown +   ts - the timestep context
40486534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
40586534af1SJed Brown 
40686534af1SJed Brown     Options Database Key:
40767b8a455SSatish Balay .    -ts_pseudo_max_dt <increment> - set pseudo max dt
40886534af1SJed Brown 
40986534af1SJed Brown     Level: advanced
41086534af1SJed Brown 
411db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
41286534af1SJed Brown @*/
4139371c9d4SSatish Balay PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) {
41486534af1SJed Brown   PetscFunctionBegin;
41586534af1SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
41686534af1SJed Brown   PetscValidLogicalCollectiveReal(ts, maxdt, 2);
417cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
41886534af1SJed Brown   PetscFunctionReturn(0);
41986534af1SJed Brown }
42086534af1SJed Brown 
42182bf6240SBarry Smith /*@
42282bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
42382bf6240SBarry Smith     is computed via the formula
42482bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
42582bf6240SBarry Smith       rather than the default update,
42682bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
42782bf6240SBarry Smith 
4283f9fe445SBarry Smith    Logically Collective on TS
42915091d37SBarry Smith 
43082bf6240SBarry Smith     Input Parameter:
43182bf6240SBarry Smith .   ts - the timestep context
43282bf6240SBarry Smith 
43382bf6240SBarry Smith     Options Database Key:
43467b8a455SSatish Balay .    -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment
43582bf6240SBarry Smith 
43615091d37SBarry Smith     Level: advanced
43715091d37SBarry Smith 
438db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
43982bf6240SBarry Smith @*/
4409371c9d4SSatish Balay PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) {
44182bf6240SBarry Smith   PetscFunctionBegin;
4420700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
443cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts));
44482bf6240SBarry Smith   PetscFunctionReturn(0);
44582bf6240SBarry Smith }
44682bf6240SBarry Smith 
447ac226902SBarry Smith /*@C
44882bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
44982bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
45082bf6240SBarry Smith 
4513f9fe445SBarry Smith    Logically Collective on TS
45215091d37SBarry Smith 
45382bf6240SBarry Smith    Input Parameters:
45415091d37SBarry Smith +  ts - timestep context
45582bf6240SBarry Smith .  dt - function to compute timestep
45615091d37SBarry Smith -  ctx - [optional] user-defined context for private data
4570298fd71SBarry Smith          required by the function (may be NULL)
45882bf6240SBarry Smith 
45915091d37SBarry Smith    Level: intermediate
460fee21e36SBarry Smith 
46182bf6240SBarry Smith    Calling sequence of func:
462a2b725a8SWilliam Gropp $  func (TS ts,PetscReal *newdt,void *ctx);
46382bf6240SBarry Smith 
464a2b725a8SWilliam Gropp +  newdt - the newly computed timestep
465a2b725a8SWilliam Gropp -  ctx - [optional] timestep context
46682bf6240SBarry Smith 
46782bf6240SBarry Smith    Notes:
46882bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
46982bf6240SBarry Smith    during the timestepping process.
4708d359177SBarry Smith    If not set then TSPseudoTimeStepDefault() is automatically used
47182bf6240SBarry Smith 
472db781477SPatrick Sanan .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
47382bf6240SBarry Smith @*/
4749371c9d4SSatish Balay PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS, PetscReal *, void *), void *ctx) {
47582bf6240SBarry Smith   PetscFunctionBegin;
4760700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
477cac4c232SBarry Smith   PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode(*)(TS, PetscReal *, void *), void *), (ts, dt, ctx));
47882bf6240SBarry Smith   PetscFunctionReturn(0);
47982bf6240SBarry Smith }
48082bf6240SBarry Smith 
48182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
48282bf6240SBarry Smith 
483ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/
4849371c9d4SSatish Balay static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) {
485be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
48682bf6240SBarry Smith 
48782bf6240SBarry Smith   PetscFunctionBegin;
48882bf6240SBarry Smith   pseudo->verify    = dt;
48982bf6240SBarry Smith   pseudo->verifyctx = ctx;
49082bf6240SBarry Smith   PetscFunctionReturn(0);
49182bf6240SBarry Smith }
49282bf6240SBarry Smith 
4939371c9d4SSatish Balay static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) {
4944bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
49582bf6240SBarry Smith 
49682bf6240SBarry Smith   PetscFunctionBegin;
49782bf6240SBarry Smith   pseudo->dt_increment = inc;
49882bf6240SBarry Smith   PetscFunctionReturn(0);
49982bf6240SBarry Smith }
50082bf6240SBarry Smith 
5019371c9d4SSatish Balay static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) {
50286534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
50386534af1SJed Brown 
50486534af1SJed Brown   PetscFunctionBegin;
50586534af1SJed Brown   pseudo->dt_max = maxdt;
50686534af1SJed Brown   PetscFunctionReturn(0);
50786534af1SJed Brown }
50886534af1SJed Brown 
5099371c9d4SSatish Balay static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) {
5104bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
51182bf6240SBarry Smith 
51282bf6240SBarry Smith   PetscFunctionBegin;
5134bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
51482bf6240SBarry Smith   PetscFunctionReturn(0);
51582bf6240SBarry Smith }
51682bf6240SBarry Smith 
5176849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/
5189371c9d4SSatish Balay static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) {
5194bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
52082bf6240SBarry Smith 
52182bf6240SBarry Smith   PetscFunctionBegin;
52282bf6240SBarry Smith   pseudo->dt    = dt;
52382bf6240SBarry Smith   pseudo->dtctx = ctx;
52482bf6240SBarry Smith   PetscFunctionReturn(0);
52582bf6240SBarry Smith }
52682bf6240SBarry Smith 
52782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
52810e6a065SJed Brown /*MC
52910e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
53082bf6240SBarry Smith 
53110e6a065SJed Brown   This method solves equations of the form
53210e6a065SJed Brown 
53310e6a065SJed Brown $    F(X,Xdot) = 0
53410e6a065SJed Brown 
53510e6a065SJed Brown   for steady state using the iteration
53610e6a065SJed Brown 
53710e6a065SJed Brown $    [G'] S = -F(X,0)
53810e6a065SJed Brown $    X += S
53910e6a065SJed Brown 
54010e6a065SJed Brown   where
54110e6a065SJed Brown 
54210e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
54310e6a065SJed Brown 
5446f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5456f2d6a7bSJed Brown   state".  See note below.
54610e6a065SJed Brown 
54710e6a065SJed Brown   Options database keys:
54810e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
5493118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
5503118ae5eSBarry Smith .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
5513118ae5eSBarry Smith -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
55210e6a065SJed Brown 
55310e6a065SJed Brown   Level: beginner
55410e6a065SJed Brown 
55510e6a065SJed Brown   References:
556606c0280SSatish Balay +  * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
557606c0280SSatish Balay -  * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
55810e6a065SJed Brown 
55910e6a065SJed Brown   Notes:
5606f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
5616f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
5626f2d6a7bSJed Brown   last step, on the first Newton iteration we have
5636f2d6a7bSJed Brown 
5646f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
5656f2d6a7bSJed Brown 
5666f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
5676f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
5686f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
56910e6a065SJed Brown 
570db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`
57110e6a065SJed Brown 
57210e6a065SJed Brown M*/
5739371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) {
5747bf11e45SBarry Smith   TS_Pseudo *pseudo;
575193ac0bcSJed Brown   SNES       snes;
57619fd82e9SBarry Smith   SNESType   stype;
5772d3f70b5SBarry Smith 
5783a40ed3dSBarry Smith   PetscFunctionBegin;
579277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
580000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
581000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
582000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
583000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
584000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
5850f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
5860f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
5872ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
588825ab935SBarry Smith   ts->usessnes            = PETSC_TRUE;
5897bf11e45SBarry Smith 
5909566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5919566063dSJacob Faibussowitsch   PetscCall(SNESGetType(snes, &stype));
5929566063dSJacob Faibussowitsch   if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY));
5932d3f70b5SBarry Smith 
594*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pseudo));
5957bf11e45SBarry Smith   ts->data = (void *)pseudo;
5962d3f70b5SBarry Smith 
597be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
598be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
59928aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6004bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
601193ac0bcSJed Brown   pseudo->fnorm                        = -1;
602be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
603be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
6043118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
6053118ae5eSBarry Smith   pseudo->fatol = 1.e-25;
6063118ae5eSBarry Smith   pseudo->frtol = 1.e-5;
6073118ae5eSBarry Smith #else
6083118ae5eSBarry Smith   pseudo->fatol = 1.e-50;
6093118ae5eSBarry Smith   pseudo->frtol = 1.e-12;
6103118ae5eSBarry Smith #endif
6119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo));
6129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo));
6139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo));
6149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo));
6159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo));
6163a40ed3dSBarry Smith   PetscFunctionReturn(0);
6172d3f70b5SBarry Smith }
6182d3f70b5SBarry Smith 
61982bf6240SBarry Smith /*@C
6208d359177SBarry Smith    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
62182bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
62228aa8177SBarry Smith 
62315091d37SBarry Smith    Collective on TS
62415091d37SBarry Smith 
62528aa8177SBarry Smith    Input Parameters:
626a2b725a8SWilliam Gropp +  ts - the timestep context
627a2b725a8SWilliam Gropp -  dtctx - unused timestep context
62828aa8177SBarry Smith 
62982bf6240SBarry Smith    Output Parameter:
63082bf6240SBarry Smith .  newdt - the timestep to use for the next step
63128aa8177SBarry Smith 
63215091d37SBarry Smith    Level: advanced
63315091d37SBarry Smith 
634db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`
63528aa8177SBarry Smith @*/
6369371c9d4SSatish Balay PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) {
63782bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo *)ts->data;
638be5899b3SLisandro Dalcin   PetscReal  inc    = pseudo->dt_increment;
63928aa8177SBarry Smith 
6403a40ed3dSBarry Smith   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(pseudo->xdot));
6429566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE));
6439566063dSJacob Faibussowitsch   PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm));
644be5899b3SLisandro Dalcin   if (pseudo->fnorm_initial < 0) {
64582bf6240SBarry Smith     /* first time through so compute initial function norm */
646cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial  = pseudo->fnorm;
647be5899b3SLisandro Dalcin     pseudo->fnorm_previous = pseudo->fnorm;
64882bf6240SBarry Smith   }
649bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step;
650bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm;
651be5899b3SLisandro Dalcin   else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm;
65286534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max);
65382bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6543a40ed3dSBarry Smith   PetscFunctionReturn(0);
65528aa8177SBarry Smith }
656