xref: /petsc/src/ts/impls/pseudo/posindep.c (revision c99c2784cfc5fd692d8686664bd1d48cf22ece0d)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4c6db04a5SJed Brown #include <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 
1887828ca2SBarry Smith   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
1987828ca2SBarry Smith   PetscReal  fnorm_previous;
2028aa8177SBarry Smith 
2187828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
22ace3abfcSBarry Smith   PetscBool  increment_dt_from_initial_dt;
237bf11e45SBarry Smith } TS_Pseudo;
242d3f70b5SBarry Smith 
252d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
262d3f70b5SBarry Smith 
274a2ae208SSatish Balay #undef __FUNCT__
284a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
297bf11e45SBarry Smith /*@
307bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
31564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
322d3f70b5SBarry Smith 
3315091d37SBarry Smith     Collective on TS
3415091d37SBarry Smith 
357bf11e45SBarry Smith     Input Parameter:
367bf11e45SBarry Smith .   ts - timestep context
377bf11e45SBarry Smith 
387bf11e45SBarry Smith     Output Parameter:
39fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
40fb4a63b6SLois Curfman McInnes 
4115091d37SBarry Smith     Level: advanced
42564e8f4eSLois Curfman McInnes 
43564e8f4eSLois Curfman McInnes     Notes:
44564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
45564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
46564e8f4eSLois Curfman McInnes 
47fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
48564e8f4eSLois Curfman McInnes 
49564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
507bf11e45SBarry Smith @*/
517087cfbeSBarry Smith PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
527bf11e45SBarry Smith {
537bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
557bf11e45SBarry Smith 
563a40ed3dSBarry Smith   PetscFunctionBegin;
57d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
587bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
59d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
617bf11e45SBarry Smith }
627bf11e45SBarry Smith 
637bf11e45SBarry Smith 
647bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
654a2ae208SSatish Balay #undef __FUNCT__
664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
677bf11e45SBarry Smith /*@C
68639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
697bf11e45SBarry Smith 
7015091d37SBarry Smith    Collective on TS
7115091d37SBarry Smith 
727bf11e45SBarry Smith    Input Parameters:
7315091d37SBarry Smith +  ts - the timestep context
747bf11e45SBarry Smith .  dtctx - unused timestep context
7515091d37SBarry Smith -  update - latest solution vector
767bf11e45SBarry Smith 
77564e8f4eSLois Curfman McInnes    Output Parameters:
7815091d37SBarry Smith +  newdt - the timestep to use for the next step
7915091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
807bf11e45SBarry Smith 
8115091d37SBarry Smith    Level: advanced
82fee21e36SBarry Smith 
83564e8f4eSLois Curfman McInnes    Note:
84564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
85564e8f4eSLois Curfman McInnes    timestep.
86564e8f4eSLois Curfman McInnes 
87564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
88564e8f4eSLois Curfman McInnes 
89564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
907bf11e45SBarry Smith @*/
917087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
927bf11e45SBarry Smith {
933a40ed3dSBarry Smith   PetscFunctionBegin;
94a7cc72afSBarry Smith   *flag = PETSC_TRUE;
953a40ed3dSBarry Smith   PetscFunctionReturn(0);
967bf11e45SBarry Smith }
977bf11e45SBarry Smith 
987bf11e45SBarry Smith 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1017bf11e45SBarry Smith /*@
102564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1037bf11e45SBarry Smith 
10415091d37SBarry Smith     Collective on TS
10515091d37SBarry Smith 
106fb4a63b6SLois Curfman McInnes     Input Parameters:
10715091d37SBarry Smith +   ts - timestep context
10815091d37SBarry Smith -   update - latest solution vector
1097bf11e45SBarry Smith 
110fb4a63b6SLois Curfman McInnes     Output Parameters:
11115091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11215091d37SBarry Smith -   flag - indicates if current timestep was ok
1137bf11e45SBarry Smith 
11415091d37SBarry Smith     Level: advanced
115fee21e36SBarry Smith 
116564e8f4eSLois Curfman McInnes     Notes:
117564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
118564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
119564e8f4eSLois Curfman McInnes 
120564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
121564e8f4eSLois Curfman McInnes 
122564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1237bf11e45SBarry Smith @*/
1247087cfbeSBarry Smith PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *flag)
1257bf11e45SBarry Smith {
1267bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
127dfbe8321SBarry Smith   PetscErrorCode ierr;
1287bf11e45SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionBegin;
130a7cc72afSBarry Smith   if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);}
1317bf11e45SBarry Smith 
1327bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1337bf11e45SBarry Smith 
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1357bf11e45SBarry Smith }
1367bf11e45SBarry Smith 
1377bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1387bf11e45SBarry Smith 
1394a2ae208SSatish Balay #undef __FUNCT__
1404a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
141193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1422d3f70b5SBarry Smith {
143277b19d0SLisandro Dalcin   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
144dfbe8321SBarry Smith   PetscErrorCode ierr;
145*c99c2784SJed Brown   PetscInt       its,lits,rej;
146ace3abfcSBarry Smith   PetscBool      ok;
147193ac0bcSJed Brown   SNESConvergedReason snesreason;
1482d3f70b5SBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
1503f2090d5SJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
151193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
152193ac0bcSJed Brown   ierr = TSPseudoComputeTimeStep(ts,&ts->next_time_step);CHKERRQ(ierr);
153193ac0bcSJed Brown   for (rej=0; ; ) {
154193ac0bcSJed Brown     ts->time_step = ts->next_time_step;
155193ac0bcSJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr);
156b850b91aSLisandro Dalcin     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1579f954729SBarry Smith     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
158c9780b6fSBarry Smith     ts->nonlinear_its += its; ts->linear_its += lits;
159193ac0bcSJed Brown     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
160193ac0bcSJed Brown     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
161193ac0bcSJed Brown     if (snesreason < 0) {
162193ac0bcSJed Brown       if (++ts->num_snes_failures >= ts->max_snes_failures) {
163193ac0bcSJed Brown         ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
164193ac0bcSJed Brown         ierr = PetscInfo2(ts,"step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
165193ac0bcSJed Brown         break;
1667bf11e45SBarry Smith       }
167193ac0bcSJed Brown     }
168193ac0bcSJed Brown     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->next_time_step,&ok);CHKERRQ(ierr);
169193ac0bcSJed Brown     if (ok) {
170193ac0bcSJed Brown       ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
171193ac0bcSJed Brown       break;
172193ac0bcSJed Brown     }
173193ac0bcSJed Brown     ts->reject++;
174193ac0bcSJed Brown     if (++rej >= ts->max_reject) {
175193ac0bcSJed Brown       ts->reason = TS_DIVERGED_STEP_REJECTED;
176193ac0bcSJed Brown     }
177193ac0bcSJed Brown   }
178193ac0bcSJed Brown   ts->ptime += ts->time_step;
1792d3f70b5SBarry Smith   ts->steps++;
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
1812d3f70b5SBarry Smith }
1822d3f70b5SBarry Smith 
1832d3f70b5SBarry Smith /*------------------------------------------------------------*/
1844a2ae208SSatish Balay #undef __FUNCT__
185277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
186277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1872d3f70b5SBarry Smith {
1887bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
189dfbe8321SBarry Smith   PetscErrorCode ierr;
1902d3f70b5SBarry Smith 
1913a40ed3dSBarry Smith   PetscFunctionBegin;
1926bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
1936bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
1946bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1962d3f70b5SBarry Smith }
1972d3f70b5SBarry Smith 
198277b19d0SLisandro Dalcin #undef __FUNCT__
199277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
200277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
201277b19d0SLisandro Dalcin {
202277b19d0SLisandro Dalcin   PetscErrorCode ierr;
203277b19d0SLisandro Dalcin 
204277b19d0SLisandro Dalcin   PetscFunctionBegin;
205277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
206277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
207277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
208277b19d0SLisandro Dalcin }
2092d3f70b5SBarry Smith 
2102d3f70b5SBarry Smith /*------------------------------------------------------------*/
2112d3f70b5SBarry Smith 
2124a2ae208SSatish Balay #undef __FUNCT__
2136f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2146f2d6a7bSJed Brown /*
2156f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2166f2d6a7bSJed Brown */
2176f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2182d3f70b5SBarry Smith {
2196f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
220193ac0bcSJed Brown   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
221193ac0bcSJed Brown   PetscScalar    *xdot;
222dfbe8321SBarry Smith   PetscErrorCode ierr;
223a7cc72afSBarry Smith   PetscInt       i,n;
2242d3f70b5SBarry Smith 
2253a40ed3dSBarry Smith   PetscFunctionBegin;
226193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
227193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2286f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2296f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2302d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2316f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2322d3f70b5SBarry Smith   }
233193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
234193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2356f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2366f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2373a40ed3dSBarry Smith   PetscFunctionReturn(0);
2382d3f70b5SBarry Smith }
2392d3f70b5SBarry Smith 
2404a2ae208SSatish Balay #undef __FUNCT__
2410f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2426f2d6a7bSJed Brown /*
2436f2d6a7bSJed Brown     The transient residual is
2446f2d6a7bSJed Brown 
2456f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2466f2d6a7bSJed Brown 
2476f2d6a7bSJed Brown     or for ODE,
2486f2d6a7bSJed Brown 
2496f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2506f2d6a7bSJed Brown 
2516f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2526f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2536f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2546f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2556f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2566f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2576f2d6a7bSJed Brown     residual, and then advances to the next time step.
2586f2d6a7bSJed Brown */
2590f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2602d3f70b5SBarry Smith {
2616f2d6a7bSJed Brown   Vec            Xdot;
262dfbe8321SBarry Smith   PetscErrorCode ierr;
2632d3f70b5SBarry Smith 
2643a40ed3dSBarry Smith   PetscFunctionBegin;
2656f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
266193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2676f2d6a7bSJed Brown   PetscFunctionReturn(0);
2686f2d6a7bSJed Brown }
2692d3f70b5SBarry Smith 
2706f2d6a7bSJed Brown #undef __FUNCT__
2710f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2726f2d6a7bSJed Brown /*
2736f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2746f2d6a7bSJed Brown 
2756f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2766f2d6a7bSJed Brown 
2776f2d6a7bSJed Brown     and for ODE:
2786f2d6a7bSJed Brown 
2796f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2806f2d6a7bSJed Brown */
2810f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2826f2d6a7bSJed Brown {
2836f2d6a7bSJed Brown   Vec            Xdot;
2846f2d6a7bSJed Brown   PetscErrorCode ierr;
2856f2d6a7bSJed Brown 
2866f2d6a7bSJed Brown   PetscFunctionBegin;
2876f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
288193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
2902d3f70b5SBarry Smith }
2912d3f70b5SBarry Smith 
2922d3f70b5SBarry Smith 
2934a2ae208SSatish Balay #undef __FUNCT__
2944a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2956849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2962d3f70b5SBarry Smith {
2977bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
298dfbe8321SBarry Smith   PetscErrorCode ierr;
2992d3f70b5SBarry Smith 
3003a40ed3dSBarry Smith   PetscFunctionBegin;
3017bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3027bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3036f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3043a40ed3dSBarry Smith   PetscFunctionReturn(0);
3052d3f70b5SBarry Smith }
3062d3f70b5SBarry Smith /*------------------------------------------------------------*/
3072d3f70b5SBarry Smith 
3084a2ae208SSatish Balay #undef __FUNCT__
309a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
310649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3112d3f70b5SBarry Smith {
3127bf11e45SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
313dfbe8321SBarry Smith   PetscErrorCode   ierr;
314649052a6SBarry Smith   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3152d3f70b5SBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
317193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
318193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
319193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
320193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
321193ac0bcSJed Brown   }
322649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
323649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
324649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3253a40ed3dSBarry Smith   PetscFunctionReturn(0);
3262d3f70b5SBarry Smith }
3272d3f70b5SBarry Smith 
3284a2ae208SSatish Balay #undef __FUNCT__
3294a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3306849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3312d3f70b5SBarry Smith {
3324bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
333dfbe8321SBarry Smith   PetscErrorCode  ierr;
334ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
335649052a6SBarry Smith   PetscViewer     viewer;
3362d3f70b5SBarry Smith 
3373a40ed3dSBarry Smith   PetscFunctionBegin;
338b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
339193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3402d3f70b5SBarry Smith     if (flg) {
341649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
342649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
34328aa8177SBarry Smith     }
34490d69ab7SBarry Smith     flg  = PETSC_FALSE;
345acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
346ca4b7087SBarry Smith     if (flg) {
347ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
348ca4b7087SBarry Smith     }
349419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
350b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3522d3f70b5SBarry Smith }
3532d3f70b5SBarry Smith 
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3566849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3572d3f70b5SBarry Smith {
3583a40ed3dSBarry Smith   PetscFunctionBegin;
3593a40ed3dSBarry Smith   PetscFunctionReturn(0);
3602d3f70b5SBarry Smith }
3612d3f70b5SBarry Smith 
36282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3634a2ae208SSatish Balay #undef __FUNCT__
3644a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
365ac226902SBarry Smith /*@C
36682bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
36782bf6240SBarry Smith    last timestep.
36882bf6240SBarry Smith 
3693f9fe445SBarry Smith    Logically Collective on TS
37015091d37SBarry Smith 
37182bf6240SBarry Smith    Input Parameters:
37215091d37SBarry Smith +  ts - timestep context
37382bf6240SBarry Smith .  dt - user-defined function to verify timestep
37415091d37SBarry Smith -  ctx - [optional] user-defined context for private data
37582bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
37682bf6240SBarry Smith 
37715091d37SBarry Smith    Level: advanced
378fee21e36SBarry Smith 
37982bf6240SBarry Smith    Calling sequence of func:
380ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
38182bf6240SBarry Smith 
38282bf6240SBarry Smith .  update - latest solution vector
38382bf6240SBarry Smith .  ctx - [optional] timestep context
38482bf6240SBarry Smith .  newdt - the timestep to use for the next step
38582bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
38682bf6240SBarry Smith 
38782bf6240SBarry Smith    Notes:
38882bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
38982bf6240SBarry Smith    during the timestepping process.
39082bf6240SBarry Smith 
39182bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
39282bf6240SBarry Smith 
39382bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
39482bf6240SBarry Smith @*/
3957087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
39682bf6240SBarry Smith {
3974ac538c5SBarry Smith   PetscErrorCode ierr;
39882bf6240SBarry Smith 
39982bf6240SBarry Smith   PetscFunctionBegin;
4000700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4014ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
40282bf6240SBarry Smith   PetscFunctionReturn(0);
40382bf6240SBarry Smith }
40482bf6240SBarry Smith 
4054a2ae208SSatish Balay #undef __FUNCT__
4064a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
40782bf6240SBarry Smith /*@
40882bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
40982bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
41082bf6240SBarry Smith 
4113f9fe445SBarry Smith    Logically Collective on TS
412fee21e36SBarry Smith 
41315091d37SBarry Smith     Input Parameters:
41415091d37SBarry Smith +   ts - the timestep context
41515091d37SBarry Smith -   inc - the scaling factor >= 1.0
41615091d37SBarry Smith 
41782bf6240SBarry Smith     Options Database Key:
41882bf6240SBarry Smith $    -ts_pseudo_increment <increment>
41982bf6240SBarry Smith 
42015091d37SBarry Smith     Level: advanced
42115091d37SBarry Smith 
42282bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42382bf6240SBarry Smith 
42482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42582bf6240SBarry Smith @*/
4267087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
42782bf6240SBarry Smith {
4284ac538c5SBarry Smith   PetscErrorCode ierr;
42982bf6240SBarry Smith 
43082bf6240SBarry Smith   PetscFunctionBegin;
4310700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
432c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4334ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
43482bf6240SBarry Smith   PetscFunctionReturn(0);
43582bf6240SBarry Smith }
43682bf6240SBarry Smith 
4374a2ae208SSatish Balay #undef __FUNCT__
4384a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
43982bf6240SBarry Smith /*@
44082bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
44182bf6240SBarry Smith     is computed via the formula
44282bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44382bf6240SBarry Smith       rather than the default update,
44482bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44582bf6240SBarry Smith 
4463f9fe445SBarry Smith    Logically Collective on TS
44715091d37SBarry Smith 
44882bf6240SBarry Smith     Input Parameter:
44982bf6240SBarry Smith .   ts - the timestep context
45082bf6240SBarry Smith 
45182bf6240SBarry Smith     Options Database Key:
45282bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45382bf6240SBarry Smith 
45415091d37SBarry Smith     Level: advanced
45515091d37SBarry Smith 
45682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
45782bf6240SBarry Smith 
45882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
45982bf6240SBarry Smith @*/
4607087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
46182bf6240SBarry Smith {
4624ac538c5SBarry Smith   PetscErrorCode ierr;
46382bf6240SBarry Smith 
46482bf6240SBarry Smith   PetscFunctionBegin;
4650700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4664ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
46782bf6240SBarry Smith   PetscFunctionReturn(0);
46882bf6240SBarry Smith }
46982bf6240SBarry Smith 
47082bf6240SBarry Smith 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
473ac226902SBarry Smith /*@C
47482bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
47582bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
47682bf6240SBarry Smith 
4773f9fe445SBarry Smith    Logically Collective on TS
47815091d37SBarry Smith 
47982bf6240SBarry Smith    Input Parameters:
48015091d37SBarry Smith +  ts - timestep context
48182bf6240SBarry Smith .  dt - function to compute timestep
48215091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48382bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48482bf6240SBarry Smith 
48515091d37SBarry Smith    Level: intermediate
486fee21e36SBarry Smith 
48782bf6240SBarry Smith    Calling sequence of func:
48887828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
48982bf6240SBarry Smith 
49082bf6240SBarry Smith .  newdt - the newly computed timestep
49182bf6240SBarry Smith .  ctx - [optional] timestep context
49282bf6240SBarry Smith 
49382bf6240SBarry Smith    Notes:
49482bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
49582bf6240SBarry Smith    during the timestepping process.
49682bf6240SBarry Smith 
49782bf6240SBarry Smith .keywords: timestep, pseudo, set
49882bf6240SBarry Smith 
49982bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
50082bf6240SBarry Smith @*/
5017087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
50282bf6240SBarry Smith {
5034ac538c5SBarry Smith   PetscErrorCode ierr;
50482bf6240SBarry Smith 
50582bf6240SBarry Smith   PetscFunctionBegin;
5060700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5074ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
50882bf6240SBarry Smith   PetscFunctionReturn(0);
50982bf6240SBarry Smith }
51082bf6240SBarry Smith 
51182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
51282bf6240SBarry Smith 
513ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
514fb2e594dSBarry Smith EXTERN_C_BEGIN
5154a2ae208SSatish Balay #undef __FUNCT__
5164a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5177087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
51882bf6240SBarry Smith {
51982bf6240SBarry Smith   TS_Pseudo *pseudo;
52082bf6240SBarry Smith 
52182bf6240SBarry Smith   PetscFunctionBegin;
52282bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
52382bf6240SBarry Smith   pseudo->verify      = dt;
52482bf6240SBarry Smith   pseudo->verifyctx   = ctx;
52582bf6240SBarry Smith   PetscFunctionReturn(0);
52682bf6240SBarry Smith }
527fb2e594dSBarry Smith EXTERN_C_END
52882bf6240SBarry Smith 
529fb2e594dSBarry Smith EXTERN_C_BEGIN
5304a2ae208SSatish Balay #undef __FUNCT__
5314a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5327087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
53382bf6240SBarry Smith {
5344bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53582bf6240SBarry Smith 
53682bf6240SBarry Smith   PetscFunctionBegin;
53782bf6240SBarry Smith   pseudo->dt_increment = inc;
53882bf6240SBarry Smith   PetscFunctionReturn(0);
53982bf6240SBarry Smith }
540fb2e594dSBarry Smith EXTERN_C_END
54182bf6240SBarry Smith 
542fb2e594dSBarry Smith EXTERN_C_BEGIN
5434a2ae208SSatish Balay #undef __FUNCT__
5444a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
5457087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
54682bf6240SBarry Smith {
5474bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54882bf6240SBarry Smith 
54982bf6240SBarry Smith   PetscFunctionBegin;
5504bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
55182bf6240SBarry Smith   PetscFunctionReturn(0);
55282bf6240SBarry Smith }
553fb2e594dSBarry Smith EXTERN_C_END
55482bf6240SBarry Smith 
5556849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
556fb2e594dSBarry Smith EXTERN_C_BEGIN
5574a2ae208SSatish Balay #undef __FUNCT__
5584a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
5597087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
56082bf6240SBarry Smith {
5614bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
56282bf6240SBarry Smith 
56382bf6240SBarry Smith   PetscFunctionBegin;
56482bf6240SBarry Smith   pseudo->dt      = dt;
56582bf6240SBarry Smith   pseudo->dtctx   = ctx;
56682bf6240SBarry Smith   PetscFunctionReturn(0);
56782bf6240SBarry Smith }
568fb2e594dSBarry Smith EXTERN_C_END
56982bf6240SBarry Smith 
57082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
57110e6a065SJed Brown /*MC
57210e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
57382bf6240SBarry Smith 
57410e6a065SJed Brown   This method solves equations of the form
57510e6a065SJed Brown 
57610e6a065SJed Brown $    F(X,Xdot) = 0
57710e6a065SJed Brown 
57810e6a065SJed Brown   for steady state using the iteration
57910e6a065SJed Brown 
58010e6a065SJed Brown $    [G'] S = -F(X,0)
58110e6a065SJed Brown $    X += S
58210e6a065SJed Brown 
58310e6a065SJed Brown   where
58410e6a065SJed Brown 
58510e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
58610e6a065SJed Brown 
5876f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5886f2d6a7bSJed Brown   state".  See note below.
58910e6a065SJed Brown 
59010e6a065SJed Brown   Options database keys:
59110e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
59210e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
59310e6a065SJed Brown 
59410e6a065SJed Brown   Level: beginner
59510e6a065SJed Brown 
59610e6a065SJed Brown   References:
59710e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
59810e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
59910e6a065SJed Brown 
60010e6a065SJed Brown   Notes:
6016f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6026f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6036f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6046f2d6a7bSJed Brown 
6056f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6066f2d6a7bSJed Brown 
6076f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6086f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6096f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
61010e6a065SJed Brown 
61110e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
61210e6a065SJed Brown 
61310e6a065SJed Brown M*/
614fb2e594dSBarry Smith EXTERN_C_BEGIN
6154a2ae208SSatish Balay #undef __FUNCT__
6164a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6177087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6182d3f70b5SBarry Smith {
6197bf11e45SBarry Smith   TS_Pseudo      *pseudo;
620dfbe8321SBarry Smith   PetscErrorCode ierr;
621193ac0bcSJed Brown   SNES           snes;
622193ac0bcSJed Brown   const SNESType stype;
6232d3f70b5SBarry Smith 
6243a40ed3dSBarry Smith   PetscFunctionBegin;
625277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
626000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
627000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6282d3f70b5SBarry Smith 
629000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
630000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
631000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6320f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6330f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6347bf11e45SBarry Smith 
635193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
636193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
637193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6382d3f70b5SBarry Smith 
63938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6407bf11e45SBarry Smith   ts->data = (void*)pseudo;
6412d3f70b5SBarry Smith 
64228aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6434bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
64428aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
645193ac0bcSJed Brown   pseudo->fnorm                        = -1;
64682bf6240SBarry Smith 
647f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
648e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6490c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
650f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
651e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6520c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
653f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
654e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6550c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6560c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6570c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6580c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6593a40ed3dSBarry Smith   PetscFunctionReturn(0);
6602d3f70b5SBarry Smith }
661fb2e594dSBarry Smith EXTERN_C_END
6622d3f70b5SBarry Smith 
6634a2ae208SSatish Balay #undef __FUNCT__
6644a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
66582bf6240SBarry Smith /*@C
66682bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
66782bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
66828aa8177SBarry Smith 
66915091d37SBarry Smith    Collective on TS
67015091d37SBarry Smith 
67128aa8177SBarry Smith    Input Parameters:
67228aa8177SBarry Smith .  ts - the timestep context
67382bf6240SBarry Smith .  dtctx - unused timestep context
67428aa8177SBarry Smith 
67582bf6240SBarry Smith    Output Parameter:
67682bf6240SBarry Smith .  newdt - the timestep to use for the next step
67728aa8177SBarry Smith 
67815091d37SBarry Smith    Level: advanced
67915091d37SBarry Smith 
68082bf6240SBarry Smith .keywords: timestep, pseudo, default
681564e8f4eSLois Curfman McInnes 
68282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
68328aa8177SBarry Smith @*/
6847087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
68528aa8177SBarry Smith {
68682bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
68787828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
688dfbe8321SBarry Smith   PetscErrorCode ierr;
68928aa8177SBarry Smith 
6903a40ed3dSBarry Smith   PetscFunctionBegin;
691bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
692214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
69382bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
69482bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
69582bf6240SBarry Smith     /* first time through so compute initial function norm */
69682bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
69782bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
69882bf6240SBarry Smith   }
69982bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
70082bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
701004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
70282bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
70382bf6240SBarry Smith   } else {
70482bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
70582bf6240SBarry Smith   }
70682bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7073a40ed3dSBarry Smith   PetscFunctionReturn(0);
70828aa8177SBarry Smith }
709