xref: /petsc/src/ts/impls/pseudo/posindep.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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 
488d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep()
497bf11e45SBarry Smith @*/
507087cfbeSBarry Smith PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
517bf11e45SBarry Smith {
527bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
537bf11e45SBarry Smith 
543a40ed3dSBarry Smith   PetscFunctionBegin;
559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0));
569566063dSJacob Faibussowitsch   PetscCall((*pseudo->dt)(ts,dt,pseudo->dtctx));
579566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0));
583a40ed3dSBarry Smith   PetscFunctionReturn(0);
597bf11e45SBarry Smith }
607bf11e45SBarry Smith 
617bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
627bf11e45SBarry Smith /*@C
638d359177SBarry Smith    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
647bf11e45SBarry Smith 
6515091d37SBarry Smith    Collective on TS
6615091d37SBarry Smith 
677bf11e45SBarry Smith    Input Parameters:
6815091d37SBarry Smith +  ts - the timestep context
697bf11e45SBarry Smith .  dtctx - unused timestep context
7015091d37SBarry Smith -  update - latest solution vector
717bf11e45SBarry Smith 
72564e8f4eSLois Curfman McInnes    Output Parameters:
7315091d37SBarry Smith +  newdt - the timestep to use for the next step
7415091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
757bf11e45SBarry Smith 
7615091d37SBarry Smith    Level: advanced
77fee21e36SBarry Smith 
78564e8f4eSLois Curfman McInnes    Note:
79564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
80564e8f4eSLois Curfman McInnes    timestep.
81564e8f4eSLois Curfman McInnes 
82564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
837bf11e45SBarry Smith @*/
848d359177SBarry Smith PetscErrorCode  TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
857bf11e45SBarry Smith {
863a40ed3dSBarry Smith   PetscFunctionBegin;
87a7cc72afSBarry Smith   *flag = PETSC_TRUE;
883a40ed3dSBarry Smith   PetscFunctionReturn(0);
897bf11e45SBarry Smith }
907bf11e45SBarry Smith 
917bf11e45SBarry Smith /*@
92564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
937bf11e45SBarry Smith 
9415091d37SBarry Smith     Collective on TS
9515091d37SBarry Smith 
96fb4a63b6SLois Curfman McInnes     Input Parameters:
9715091d37SBarry Smith +   ts - timestep context
9815091d37SBarry Smith -   update - latest solution vector
997bf11e45SBarry Smith 
100fb4a63b6SLois Curfman McInnes     Output Parameters:
10115091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
10215091d37SBarry Smith -   flag - indicates if current timestep was ok
1037bf11e45SBarry Smith 
10415091d37SBarry Smith     Level: advanced
105fee21e36SBarry Smith 
106564e8f4eSLois Curfman McInnes     Notes:
107564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
108564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
109564e8f4eSLois Curfman McInnes 
1108d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
1117bf11e45SBarry Smith @*/
1127087cfbeSBarry Smith PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag)
1137bf11e45SBarry Smith {
1147bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
1157bf11e45SBarry Smith 
1163a40ed3dSBarry Smith   PetscFunctionBegin;
117cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
118be5899b3SLisandro Dalcin   if (pseudo->verify) {
1199566063dSJacob Faibussowitsch     PetscCall((*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag));
120cb9d8021SPierre Barbier de Reuille   }
1213a40ed3dSBarry Smith   PetscFunctionReturn(0);
1227bf11e45SBarry Smith }
1237bf11e45SBarry Smith 
1247bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1257bf11e45SBarry Smith 
126193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1272d3f70b5SBarry Smith {
128277b19d0SLisandro Dalcin   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
129be5899b3SLisandro Dalcin   PetscInt            nits,lits,reject;
130cdbf8f93SLisandro Dalcin   PetscBool           stepok;
131be5899b3SLisandro Dalcin   PetscReal           next_time_step = ts->time_step;
1322d3f70b5SBarry Smith 
1333a40ed3dSBarry Smith   PetscFunctionBegin;
134bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
1359566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,pseudo->update));
1369566063dSJacob Faibussowitsch   PetscCall(TSPseudoComputeTimeStep(ts,&next_time_step));
137cdbf8f93SLisandro Dalcin   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
138cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
1399566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,ts->ptime+ts->time_step));
1409566063dSJacob Faibussowitsch     PetscCall(SNESSolve(ts->snes,NULL,pseudo->update));
1419566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(ts->snes,&nits));
1429566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits));
143be5899b3SLisandro Dalcin     ts->snes_its += nits; ts->ksp_its += lits;
1449566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update)));
1459566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok));
146be5899b3SLisandro Dalcin     if (!stepok) {next_time_step = ts->time_step; continue;}
147193ac0bcSJed Brown     pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
1489566063dSJacob Faibussowitsch     PetscCall(TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok));
149cdbf8f93SLisandro Dalcin     if (stepok) break;
150cdbf8f93SLisandro Dalcin   }
151be5899b3SLisandro Dalcin   if (reject >= ts->max_reject) {
152be5899b3SLisandro Dalcin     ts->reason = TS_DIVERGED_STEP_REJECTED;
1539566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject));
154cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1557bf11e45SBarry Smith   }
156be5899b3SLisandro Dalcin 
1579566063dSJacob Faibussowitsch   PetscCall(VecCopy(pseudo->update,ts->vec_sol));
158be5899b3SLisandro Dalcin   ts->ptime += ts->time_step;
159be5899b3SLisandro Dalcin   ts->time_step = next_time_step;
160be5899b3SLisandro Dalcin 
1613118ae5eSBarry Smith   if (pseudo->fnorm < 0) {
1629566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
1639566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
1649566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
1653118ae5eSBarry Smith   }
1663118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1673118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
1689566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ts,"Step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol));
1693118ae5eSBarry Smith     PetscFunctionReturn(0);
1703118ae5eSBarry Smith   }
1713118ae5eSBarry Smith   if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) {
1723118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
1739566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ts,"Step=%D, converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,pseudo->fnorm,pseudo->fnorm_initial,pseudo->fatol));
1743118ae5eSBarry Smith     PetscFunctionReturn(0);
1753118ae5eSBarry Smith   }
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1772d3f70b5SBarry Smith }
1782d3f70b5SBarry Smith 
1792d3f70b5SBarry Smith /*------------------------------------------------------------*/
180277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1812d3f70b5SBarry Smith {
1827bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
1832d3f70b5SBarry Smith 
1843a40ed3dSBarry Smith   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->update));
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->func));
1879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pseudo->xdot));
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1892d3f70b5SBarry Smith }
1902d3f70b5SBarry Smith 
191277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
192277b19d0SLisandro Dalcin {
193277b19d0SLisandro Dalcin   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(TSReset_Pseudo(ts));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
1969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL));
1979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL));
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL));
1999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL));
201277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
202277b19d0SLisandro Dalcin }
2032d3f70b5SBarry Smith 
2042d3f70b5SBarry Smith /*------------------------------------------------------------*/
2052d3f70b5SBarry Smith 
2066f2d6a7bSJed Brown /*
2076f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2086f2d6a7bSJed Brown */
2096f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2102d3f70b5SBarry Smith {
2116f2d6a7bSJed Brown   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
212193ac0bcSJed Brown   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
213193ac0bcSJed Brown   PetscScalar       *xdot;
214a7cc72afSBarry Smith   PetscInt          i,n;
2152d3f70b5SBarry Smith 
2163a40ed3dSBarry Smith   PetscFunctionBegin;
217aab5bcd8SJed Brown   *Xdot = NULL;
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ts->vec_sol,&xn));
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&xnp1));
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pseudo->xdot,&xdot));
2219566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
222bbd56ea5SKarl Rupp   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ts->vec_sol,&xn));
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&xnp1));
2259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pseudo->xdot,&xdot));
2266f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
2282d3f70b5SBarry Smith }
2292d3f70b5SBarry Smith 
2306f2d6a7bSJed Brown /*
2316f2d6a7bSJed Brown     The transient residual is
2326f2d6a7bSJed Brown 
2336f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2346f2d6a7bSJed Brown 
2356f2d6a7bSJed Brown     or for ODE,
2366f2d6a7bSJed Brown 
2376f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2386f2d6a7bSJed Brown 
2396f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2406f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2416f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2426f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2436f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2446f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2456f2d6a7bSJed Brown     residual, and then advances to the next time step.
2466f2d6a7bSJed Brown */
2470f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2482d3f70b5SBarry Smith {
2496f2d6a7bSJed Brown   Vec            Xdot;
2502d3f70b5SBarry Smith 
2513a40ed3dSBarry Smith   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts,X,&Xdot));
2539566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE));
2546f2d6a7bSJed Brown   PetscFunctionReturn(0);
2556f2d6a7bSJed Brown }
2562d3f70b5SBarry Smith 
2576f2d6a7bSJed Brown /*
2586f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2596f2d6a7bSJed Brown 
2606f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2616f2d6a7bSJed Brown 
2626f2d6a7bSJed Brown     and for ODE:
2636f2d6a7bSJed Brown 
2646f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2656f2d6a7bSJed Brown */
266d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
2676f2d6a7bSJed Brown {
2686f2d6a7bSJed Brown   Vec            Xdot;
2696f2d6a7bSJed Brown 
2706f2d6a7bSJed Brown   PetscFunctionBegin;
2719566063dSJacob Faibussowitsch   PetscCall(TSPseudoGetXdot(ts,X,&Xdot));
2729566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE));
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2742d3f70b5SBarry Smith }
2752d3f70b5SBarry Smith 
2766849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2772d3f70b5SBarry Smith {
2787bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
2792d3f70b5SBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
2819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->update));
2829566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->func));
2839566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->xdot));
2843a40ed3dSBarry Smith   PetscFunctionReturn(0);
2852d3f70b5SBarry Smith }
2862d3f70b5SBarry Smith /*------------------------------------------------------------*/
2872d3f70b5SBarry Smith 
288560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
2892d3f70b5SBarry Smith {
2907bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
291ce94432eSBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
2922d3f70b5SBarry Smith 
2933a40ed3dSBarry Smith   PetscFunctionBegin;
294193ac0bcSJed Brown   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
2959566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(pseudo->xdot));
2969566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
2979566063dSJacob Faibussowitsch     PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
298193ac0bcSJed Brown   }
2999566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel));
3009566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm));
3019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel));
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
3032d3f70b5SBarry Smith }
3042d3f70b5SBarry Smith 
3054416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts)
3062d3f70b5SBarry Smith {
3074bbc92c1SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
308ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
309649052a6SBarry Smith   PetscViewer    viewer;
3102d3f70b5SBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options"));
3139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL));
3142d3f70b5SBarry Smith   if (flg) {
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer));
3169566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
31728aa8177SBarry Smith   }
318be5899b3SLisandro Dalcin   flg  = pseudo->increment_dt_from_initial_dt;
3199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL));
320be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
3219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL));
3229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL));
3239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL));
3249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL));
3259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
3272d3f70b5SBarry Smith }
3282d3f70b5SBarry Smith 
3296849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3302d3f70b5SBarry Smith {
3313118ae5eSBarry Smith   PetscBool      isascii;
332d52bd9f3SBarry Smith 
3333a40ed3dSBarry Smith   PetscFunctionBegin;
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
3353118ae5eSBarry Smith   if (isascii) {
3363118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  frtol - relative tolerance in function value %g\n",(double)pseudo->frtol));
3389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol));
3399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_initial - initial timestep %g\n",(double)pseudo->dt_initial));
3409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment));
3419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_max - maximum time %g\n",(double)pseudo->dt_max));
3423118ae5eSBarry Smith   }
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
3442d3f70b5SBarry Smith }
3452d3f70b5SBarry Smith 
34682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
347ac226902SBarry Smith /*@C
34882bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
34982bf6240SBarry Smith    last timestep.
35082bf6240SBarry Smith 
3513f9fe445SBarry Smith    Logically Collective on TS
35215091d37SBarry Smith 
35382bf6240SBarry Smith    Input Parameters:
35415091d37SBarry Smith +  ts - timestep context
35582bf6240SBarry Smith .  dt - user-defined function to verify timestep
35615091d37SBarry Smith -  ctx - [optional] user-defined context for private data
3570298fd71SBarry Smith          for the timestep verification routine (may be NULL)
35882bf6240SBarry Smith 
35915091d37SBarry Smith    Level: advanced
360fee21e36SBarry Smith 
36182bf6240SBarry Smith    Calling sequence of func:
362a2b725a8SWilliam Gropp $  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
36382bf6240SBarry Smith 
364a2b725a8SWilliam Gropp +  update - latest solution vector
36582bf6240SBarry Smith .  ctx - [optional] timestep context
36682bf6240SBarry Smith .  newdt - the timestep to use for the next step
367a2b725a8SWilliam Gropp -  flag - flag indicating whether the last time step was acceptable
36882bf6240SBarry Smith 
36982bf6240SBarry Smith    Notes:
37082bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
37182bf6240SBarry Smith    during the timestepping process.
37282bf6240SBarry Smith 
3738d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
37482bf6240SBarry Smith @*/
3757087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
37682bf6240SBarry Smith {
37782bf6240SBarry Smith   PetscFunctionBegin;
3780700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
379*cac4c232SBarry Smith   PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));
38082bf6240SBarry Smith   PetscFunctionReturn(0);
38182bf6240SBarry Smith }
38282bf6240SBarry Smith 
38382bf6240SBarry Smith /*@
38482bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
3858d359177SBarry Smith     dt when using the TSPseudoTimeStepDefault() routine.
38682bf6240SBarry Smith 
3873f9fe445SBarry Smith    Logically Collective on TS
388fee21e36SBarry Smith 
38915091d37SBarry Smith     Input Parameters:
39015091d37SBarry Smith +   ts - the timestep context
39115091d37SBarry Smith -   inc - the scaling factor >= 1.0
39215091d37SBarry Smith 
39382bf6240SBarry Smith     Options Database Key:
39467b8a455SSatish Balay .    -ts_pseudo_increment <increment> - set pseudo increment
39582bf6240SBarry Smith 
39615091d37SBarry Smith     Level: advanced
39715091d37SBarry Smith 
3988d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
39982bf6240SBarry Smith @*/
4007087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
40182bf6240SBarry Smith {
40282bf6240SBarry Smith   PetscFunctionBegin;
4030700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
404c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
405*cac4c232SBarry Smith   PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));
40682bf6240SBarry Smith   PetscFunctionReturn(0);
40782bf6240SBarry Smith }
40882bf6240SBarry Smith 
40986534af1SJed Brown /*@
41086534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
4118d359177SBarry Smith     when using the TSPseudoTimeStepDefault() routine.
41286534af1SJed Brown 
41386534af1SJed Brown    Logically Collective on TS
41486534af1SJed Brown 
41586534af1SJed Brown     Input Parameters:
41686534af1SJed Brown +   ts - the timestep context
41786534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
41886534af1SJed Brown 
41986534af1SJed Brown     Options Database Key:
42067b8a455SSatish Balay .    -ts_pseudo_max_dt <increment> - set pseudo max dt
42186534af1SJed Brown 
42286534af1SJed Brown     Level: advanced
42386534af1SJed Brown 
4248d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
42586534af1SJed Brown @*/
42686534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
42786534af1SJed Brown {
42886534af1SJed Brown   PetscFunctionBegin;
42986534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
43086534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
431*cac4c232SBarry Smith   PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
43286534af1SJed Brown   PetscFunctionReturn(0);
43386534af1SJed Brown }
43486534af1SJed Brown 
43582bf6240SBarry Smith /*@
43682bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
43782bf6240SBarry Smith     is computed via the formula
43882bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
43982bf6240SBarry Smith       rather than the default update,
44082bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44182bf6240SBarry Smith 
4423f9fe445SBarry Smith    Logically Collective on TS
44315091d37SBarry Smith 
44482bf6240SBarry Smith     Input Parameter:
44582bf6240SBarry Smith .   ts - the timestep context
44682bf6240SBarry Smith 
44782bf6240SBarry Smith     Options Database Key:
44867b8a455SSatish Balay .    -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment
44982bf6240SBarry Smith 
45015091d37SBarry Smith     Level: advanced
45115091d37SBarry Smith 
4528d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
45382bf6240SBarry Smith @*/
4547087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
45582bf6240SBarry Smith {
45682bf6240SBarry Smith   PetscFunctionBegin;
4570700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
458*cac4c232SBarry Smith   PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));
45982bf6240SBarry Smith   PetscFunctionReturn(0);
46082bf6240SBarry Smith }
46182bf6240SBarry Smith 
462ac226902SBarry Smith /*@C
46382bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
46482bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
46582bf6240SBarry Smith 
4663f9fe445SBarry Smith    Logically Collective on TS
46715091d37SBarry Smith 
46882bf6240SBarry Smith    Input Parameters:
46915091d37SBarry Smith +  ts - timestep context
47082bf6240SBarry Smith .  dt - function to compute timestep
47115091d37SBarry Smith -  ctx - [optional] user-defined context for private data
4720298fd71SBarry Smith          required by the function (may be NULL)
47382bf6240SBarry Smith 
47415091d37SBarry Smith    Level: intermediate
475fee21e36SBarry Smith 
47682bf6240SBarry Smith    Calling sequence of func:
477a2b725a8SWilliam Gropp $  func (TS ts,PetscReal *newdt,void *ctx);
47882bf6240SBarry Smith 
479a2b725a8SWilliam Gropp +  newdt - the newly computed timestep
480a2b725a8SWilliam Gropp -  ctx - [optional] timestep context
48182bf6240SBarry Smith 
48282bf6240SBarry Smith    Notes:
48382bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
48482bf6240SBarry Smith    during the timestepping process.
4858d359177SBarry Smith    If not set then TSPseudoTimeStepDefault() is automatically used
48682bf6240SBarry Smith 
4878d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
48882bf6240SBarry Smith @*/
4897087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
49082bf6240SBarry Smith {
49182bf6240SBarry Smith   PetscFunctionBegin;
4920700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
493*cac4c232SBarry Smith   PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));
49482bf6240SBarry Smith   PetscFunctionReturn(0);
49582bf6240SBarry Smith }
49682bf6240SBarry Smith 
49782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
49882bf6240SBarry Smith 
499ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
500560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
50182bf6240SBarry Smith {
502be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
50382bf6240SBarry Smith 
50482bf6240SBarry Smith   PetscFunctionBegin;
50582bf6240SBarry Smith   pseudo->verify    = dt;
50682bf6240SBarry Smith   pseudo->verifyctx = ctx;
50782bf6240SBarry Smith   PetscFunctionReturn(0);
50882bf6240SBarry Smith }
50982bf6240SBarry Smith 
510560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
51182bf6240SBarry Smith {
5124bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
51382bf6240SBarry Smith 
51482bf6240SBarry Smith   PetscFunctionBegin;
51582bf6240SBarry Smith   pseudo->dt_increment = inc;
51682bf6240SBarry Smith   PetscFunctionReturn(0);
51782bf6240SBarry Smith }
51882bf6240SBarry Smith 
519560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
52086534af1SJed Brown {
52186534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
52286534af1SJed Brown 
52386534af1SJed Brown   PetscFunctionBegin;
52486534af1SJed Brown   pseudo->dt_max = maxdt;
52586534af1SJed Brown   PetscFunctionReturn(0);
52686534af1SJed Brown }
52786534af1SJed Brown 
528560360afSLisandro Dalcin static PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
52982bf6240SBarry Smith {
5304bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53182bf6240SBarry Smith 
53282bf6240SBarry Smith   PetscFunctionBegin;
5334bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
53482bf6240SBarry Smith   PetscFunctionReturn(0);
53582bf6240SBarry Smith }
53682bf6240SBarry Smith 
5376849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
538560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
53982bf6240SBarry Smith {
5404bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54182bf6240SBarry Smith 
54282bf6240SBarry Smith   PetscFunctionBegin;
54382bf6240SBarry Smith   pseudo->dt    = dt;
54482bf6240SBarry Smith   pseudo->dtctx = ctx;
54582bf6240SBarry Smith   PetscFunctionReturn(0);
54682bf6240SBarry Smith }
54782bf6240SBarry Smith 
54882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
54910e6a065SJed Brown /*MC
55010e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
55182bf6240SBarry Smith 
55210e6a065SJed Brown   This method solves equations of the form
55310e6a065SJed Brown 
55410e6a065SJed Brown $    F(X,Xdot) = 0
55510e6a065SJed Brown 
55610e6a065SJed Brown   for steady state using the iteration
55710e6a065SJed Brown 
55810e6a065SJed Brown $    [G'] S = -F(X,0)
55910e6a065SJed Brown $    X += S
56010e6a065SJed Brown 
56110e6a065SJed Brown   where
56210e6a065SJed Brown 
56310e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
56410e6a065SJed Brown 
5656f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5666f2d6a7bSJed Brown   state".  See note below.
56710e6a065SJed Brown 
56810e6a065SJed Brown   Options database keys:
56910e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
5703118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
5713118ae5eSBarry Smith .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
5723118ae5eSBarry Smith -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
57310e6a065SJed Brown 
57410e6a065SJed Brown   Level: beginner
57510e6a065SJed Brown 
57610e6a065SJed Brown   References:
577606c0280SSatish Balay +  * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
578606c0280SSatish Balay -  * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
57910e6a065SJed Brown 
58010e6a065SJed Brown   Notes:
5816f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
5826f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
5836f2d6a7bSJed Brown   last step, on the first Newton iteration we have
5846f2d6a7bSJed Brown 
5856f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
5866f2d6a7bSJed Brown 
5876f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
5886f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
5896f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
59010e6a065SJed Brown 
59110e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
59210e6a065SJed Brown 
59310e6a065SJed Brown M*/
5948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
5952d3f70b5SBarry Smith {
5967bf11e45SBarry Smith   TS_Pseudo      *pseudo;
597193ac0bcSJed Brown   SNES           snes;
59819fd82e9SBarry Smith   SNESType       stype;
5992d3f70b5SBarry Smith 
6003a40ed3dSBarry Smith   PetscFunctionBegin;
601277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
602000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
603000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
604000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
605000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
606000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
6070f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
6080f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
6092ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
610825ab935SBarry Smith   ts->usessnes            = PETSC_TRUE;
6117bf11e45SBarry Smith 
6129566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
6139566063dSJacob Faibussowitsch   PetscCall(SNESGetType(snes,&stype));
6149566063dSJacob Faibussowitsch   if (!stype) PetscCall(SNESSetType(snes,SNESKSPONLY));
6152d3f70b5SBarry Smith 
6169566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&pseudo));
6177bf11e45SBarry Smith   ts->data = (void*)pseudo;
6182d3f70b5SBarry Smith 
619be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
620be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
62128aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6224bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
623193ac0bcSJed Brown   pseudo->fnorm                        = -1;
624be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
625be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
6263118ae5eSBarry Smith  #if defined(PETSC_USE_REAL_SINGLE)
6273118ae5eSBarry Smith   pseudo->fatol                        = 1.e-25;
6283118ae5eSBarry Smith   pseudo->frtol                        = 1.e-5;
6293118ae5eSBarry Smith #else
6303118ae5eSBarry Smith   pseudo->fatol                        = 1.e-50;
6313118ae5eSBarry Smith   pseudo->frtol                        = 1.e-12;
6323118ae5eSBarry Smith #endif
6339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo));
6349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo));
6359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo));
6369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo));
6379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo));
6383a40ed3dSBarry Smith   PetscFunctionReturn(0);
6392d3f70b5SBarry Smith }
6402d3f70b5SBarry Smith 
64182bf6240SBarry Smith /*@C
6428d359177SBarry Smith    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
64382bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
64428aa8177SBarry Smith 
64515091d37SBarry Smith    Collective on TS
64615091d37SBarry Smith 
64728aa8177SBarry Smith    Input Parameters:
648a2b725a8SWilliam Gropp +  ts - the timestep context
649a2b725a8SWilliam Gropp -  dtctx - unused timestep context
65028aa8177SBarry Smith 
65182bf6240SBarry Smith    Output Parameter:
65282bf6240SBarry Smith .  newdt - the timestep to use for the next step
65328aa8177SBarry Smith 
65415091d37SBarry Smith    Level: advanced
65515091d37SBarry Smith 
65682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
65728aa8177SBarry Smith @*/
6588d359177SBarry Smith PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
65928aa8177SBarry Smith {
66082bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
661be5899b3SLisandro Dalcin   PetscReal      inc = pseudo->dt_increment;
66228aa8177SBarry Smith 
6633a40ed3dSBarry Smith   PetscFunctionBegin;
6649566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(pseudo->xdot));
6659566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
6669566063dSJacob Faibussowitsch   PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
667be5899b3SLisandro Dalcin   if (pseudo->fnorm_initial < 0) {
66882bf6240SBarry Smith     /* first time through so compute initial function norm */
669cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial  = pseudo->fnorm;
670be5899b3SLisandro Dalcin     pseudo->fnorm_previous = pseudo->fnorm;
67182bf6240SBarry Smith   }
672bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
673bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
674be5899b3SLisandro Dalcin   else                                           *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm;
67586534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
67682bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6773a40ed3dSBarry Smith   PetscFunctionReturn(0);
67828aa8177SBarry Smith }
679