xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 214bc6a2d8e2ea6a1f725944c2fd032dc3e2c3e7)
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"
141a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime)
1422d3f70b5SBarry Smith {
143277b19d0SLisandro Dalcin   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
144277b19d0SLisandro Dalcin   Vec            sol = ts->vec_sol, update = pseudo->update;
145dfbe8321SBarry Smith   PetscErrorCode ierr;
146186e87acSLisandro Dalcin   PetscInt       i,its,lits;
147ace3abfcSBarry Smith   PetscBool      ok;
14887828ca2SBarry Smith   PetscReal      current_time_step;
1492d3f70b5SBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1512d3f70b5SBarry Smith   *steps = -ts->steps;
152186e87acSLisandro Dalcin   *ptime =  ts->ptime;
1532d3f70b5SBarry Smith 
154277b19d0SLisandro Dalcin   ierr = VecCopy(sol,update);CHKERRQ(ierr);
155186e87acSLisandro Dalcin   for (i=0; i<ts->max_steps && ts->ptime < ts->max_time; i++) {
156186e87acSLisandro Dalcin 
1577bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
158e6e5fe25SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1597bf11e45SBarry Smith     current_time_step = ts->time_step;
1603f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
161e6e5fe25SBarry Smith     while (PETSC_TRUE) {
1627bf11e45SBarry Smith       ts->ptime  += current_time_step;
163277b19d0SLisandro Dalcin       ierr = SNESSolve(ts->snes,PETSC_NULL,update);CHKERRQ(ierr);
164b850b91aSLisandro Dalcin       ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1659f954729SBarry Smith       ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
166c9780b6fSBarry Smith       ts->nonlinear_its += its; ts->linear_its += lits;
167277b19d0SLisandro Dalcin       ierr = TSPseudoVerifyTimeStep(ts,update,&ts->time_step,&ok);CHKERRQ(ierr);
1687bf11e45SBarry Smith       if (ok) break;
1697bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1707bf11e45SBarry Smith       current_time_step = ts->time_step;
1717bf11e45SBarry Smith     }
172277b19d0SLisandro Dalcin     ierr = VecCopy(update,sol);CHKERRQ(ierr);
1732d3f70b5SBarry Smith     ts->steps++;
1743f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
1752d3f70b5SBarry Smith   }
176bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
177*214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
178e6e5fe25SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
179e6e5fe25SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1802d3f70b5SBarry Smith 
1812d3f70b5SBarry Smith   *steps += ts->steps;
182142b95e3SSatish Balay   *ptime  = ts->ptime;
1833a40ed3dSBarry Smith   PetscFunctionReturn(0);
1842d3f70b5SBarry Smith }
1852d3f70b5SBarry Smith 
1862d3f70b5SBarry Smith /*------------------------------------------------------------*/
1874a2ae208SSatish Balay #undef __FUNCT__
188277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
189277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1902d3f70b5SBarry Smith {
1917bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
192dfbe8321SBarry Smith   PetscErrorCode ierr;
1932d3f70b5SBarry Smith 
1943a40ed3dSBarry Smith   PetscFunctionBegin;
1956bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
1966bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
1976bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1992d3f70b5SBarry Smith }
2002d3f70b5SBarry Smith 
201277b19d0SLisandro Dalcin #undef __FUNCT__
202277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
203277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
204277b19d0SLisandro Dalcin {
205277b19d0SLisandro Dalcin   PetscErrorCode ierr;
206277b19d0SLisandro Dalcin 
207277b19d0SLisandro Dalcin   PetscFunctionBegin;
208277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
209277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
210277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
211277b19d0SLisandro Dalcin }
2122d3f70b5SBarry Smith 
2132d3f70b5SBarry Smith /*------------------------------------------------------------*/
2142d3f70b5SBarry Smith 
2154a2ae208SSatish Balay #undef __FUNCT__
2166f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2176f2d6a7bSJed Brown /*
2186f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2196f2d6a7bSJed Brown */
2206f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2212d3f70b5SBarry Smith {
2226f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
2236f2d6a7bSJed Brown   PetscScalar    mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot;
224dfbe8321SBarry Smith   PetscErrorCode ierr;
225a7cc72afSBarry Smith   PetscInt       i,n;
2262d3f70b5SBarry Smith 
2273a40ed3dSBarry Smith   PetscFunctionBegin;
2286f2d6a7bSJed Brown   ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2296f2d6a7bSJed Brown   ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr);
2306f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2316f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2322d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2336f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2342d3f70b5SBarry Smith   }
2356f2d6a7bSJed Brown   ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2366f2d6a7bSJed Brown   ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr);
2376f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2386f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
2402d3f70b5SBarry Smith }
2412d3f70b5SBarry Smith 
2424a2ae208SSatish Balay #undef __FUNCT__
2430f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2446f2d6a7bSJed Brown /*
2456f2d6a7bSJed Brown     The transient residual is
2466f2d6a7bSJed Brown 
2476f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2486f2d6a7bSJed Brown 
2496f2d6a7bSJed Brown     or for ODE,
2506f2d6a7bSJed Brown 
2516f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2526f2d6a7bSJed Brown 
2536f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2546f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2556f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2566f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2576f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2586f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2596f2d6a7bSJed Brown     residual, and then advances to the next time step.
2606f2d6a7bSJed Brown */
2610f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2622d3f70b5SBarry Smith {
2636f2d6a7bSJed Brown   Vec            Xdot;
264dfbe8321SBarry Smith   PetscErrorCode ierr;
2652d3f70b5SBarry Smith 
2663a40ed3dSBarry Smith   PetscFunctionBegin;
2676f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
268*214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2696f2d6a7bSJed Brown   PetscFunctionReturn(0);
2706f2d6a7bSJed Brown }
2712d3f70b5SBarry Smith 
2726f2d6a7bSJed Brown #undef __FUNCT__
2730f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2746f2d6a7bSJed Brown /*
2756f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2766f2d6a7bSJed Brown 
2776f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2786f2d6a7bSJed Brown 
2796f2d6a7bSJed Brown     and for ODE:
2806f2d6a7bSJed Brown 
2816f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2826f2d6a7bSJed Brown */
2830f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2846f2d6a7bSJed Brown {
2856f2d6a7bSJed Brown   Vec            Xdot;
2866f2d6a7bSJed Brown   PetscErrorCode ierr;
2876f2d6a7bSJed Brown 
2886f2d6a7bSJed Brown   PetscFunctionBegin;
2896f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
290*214bc6a2SJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
2913a40ed3dSBarry Smith   PetscFunctionReturn(0);
2922d3f70b5SBarry Smith }
2932d3f70b5SBarry Smith 
2942d3f70b5SBarry Smith 
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2976849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2982d3f70b5SBarry Smith {
2997bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
300dfbe8321SBarry Smith   PetscErrorCode ierr;
3012d3f70b5SBarry Smith 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3037bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3047bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3056f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
3072d3f70b5SBarry Smith }
3082d3f70b5SBarry Smith /*------------------------------------------------------------*/
3092d3f70b5SBarry Smith 
3104a2ae208SSatish Balay #undef __FUNCT__
311a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
312649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3132d3f70b5SBarry Smith {
3147bf11e45SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
315dfbe8321SBarry Smith   PetscErrorCode   ierr;
316649052a6SBarry Smith   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3172d3f70b5SBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
319649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
320649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
321649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
3232d3f70b5SBarry Smith }
3242d3f70b5SBarry Smith 
3254a2ae208SSatish Balay #undef __FUNCT__
3264a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3276849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3282d3f70b5SBarry Smith {
3294bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
330dfbe8321SBarry Smith   PetscErrorCode  ierr;
331ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
332649052a6SBarry Smith   PetscViewer     viewer;
3332d3f70b5SBarry Smith 
3343a40ed3dSBarry Smith   PetscFunctionBegin;
335b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
336acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3372d3f70b5SBarry Smith     if (flg) {
338649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
339649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
34028aa8177SBarry Smith     }
34190d69ab7SBarry Smith     flg  = PETSC_FALSE;
342acfcf0e5SJed 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);
343ca4b7087SBarry Smith     if (flg) {
344ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
345ca4b7087SBarry Smith     }
346419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
347b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3483a40ed3dSBarry Smith   PetscFunctionReturn(0);
3492d3f70b5SBarry Smith }
3502d3f70b5SBarry Smith 
3514a2ae208SSatish Balay #undef __FUNCT__
3524a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3536849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3542d3f70b5SBarry Smith {
3553a40ed3dSBarry Smith   PetscFunctionBegin;
3563a40ed3dSBarry Smith   PetscFunctionReturn(0);
3572d3f70b5SBarry Smith }
3582d3f70b5SBarry Smith 
35982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3604a2ae208SSatish Balay #undef __FUNCT__
3614a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
362ac226902SBarry Smith /*@C
36382bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
36482bf6240SBarry Smith    last timestep.
36582bf6240SBarry Smith 
3663f9fe445SBarry Smith    Logically Collective on TS
36715091d37SBarry Smith 
36882bf6240SBarry Smith    Input Parameters:
36915091d37SBarry Smith +  ts - timestep context
37082bf6240SBarry Smith .  dt - user-defined function to verify timestep
37115091d37SBarry Smith -  ctx - [optional] user-defined context for private data
37282bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
37382bf6240SBarry Smith 
37415091d37SBarry Smith    Level: advanced
375fee21e36SBarry Smith 
37682bf6240SBarry Smith    Calling sequence of func:
377ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
37882bf6240SBarry Smith 
37982bf6240SBarry Smith .  update - latest solution vector
38082bf6240SBarry Smith .  ctx - [optional] timestep context
38182bf6240SBarry Smith .  newdt - the timestep to use for the next step
38282bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
38382bf6240SBarry Smith 
38482bf6240SBarry Smith    Notes:
38582bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
38682bf6240SBarry Smith    during the timestepping process.
38782bf6240SBarry Smith 
38882bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
38982bf6240SBarry Smith 
39082bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
39182bf6240SBarry Smith @*/
3927087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
39382bf6240SBarry Smith {
3944ac538c5SBarry Smith   PetscErrorCode ierr;
39582bf6240SBarry Smith 
39682bf6240SBarry Smith   PetscFunctionBegin;
3970700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3984ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
39982bf6240SBarry Smith   PetscFunctionReturn(0);
40082bf6240SBarry Smith }
40182bf6240SBarry Smith 
4024a2ae208SSatish Balay #undef __FUNCT__
4034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
40482bf6240SBarry Smith /*@
40582bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
40682bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
40782bf6240SBarry Smith 
4083f9fe445SBarry Smith    Logically Collective on TS
409fee21e36SBarry Smith 
41015091d37SBarry Smith     Input Parameters:
41115091d37SBarry Smith +   ts - the timestep context
41215091d37SBarry Smith -   inc - the scaling factor >= 1.0
41315091d37SBarry Smith 
41482bf6240SBarry Smith     Options Database Key:
41582bf6240SBarry Smith $    -ts_pseudo_increment <increment>
41682bf6240SBarry Smith 
41715091d37SBarry Smith     Level: advanced
41815091d37SBarry Smith 
41982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42082bf6240SBarry Smith 
42182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42282bf6240SBarry Smith @*/
4237087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
42482bf6240SBarry Smith {
4254ac538c5SBarry Smith   PetscErrorCode ierr;
42682bf6240SBarry Smith 
42782bf6240SBarry Smith   PetscFunctionBegin;
4280700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
429c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4304ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
43182bf6240SBarry Smith   PetscFunctionReturn(0);
43282bf6240SBarry Smith }
43382bf6240SBarry Smith 
4344a2ae208SSatish Balay #undef __FUNCT__
4354a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
43682bf6240SBarry Smith /*@
43782bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
43882bf6240SBarry Smith     is computed via the formula
43982bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44082bf6240SBarry Smith       rather than the default update,
44182bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44282bf6240SBarry Smith 
4433f9fe445SBarry Smith    Logically Collective on TS
44415091d37SBarry Smith 
44582bf6240SBarry Smith     Input Parameter:
44682bf6240SBarry Smith .   ts - the timestep context
44782bf6240SBarry Smith 
44882bf6240SBarry Smith     Options Database Key:
44982bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45082bf6240SBarry Smith 
45115091d37SBarry Smith     Level: advanced
45215091d37SBarry Smith 
45382bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
45482bf6240SBarry Smith 
45582bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
45682bf6240SBarry Smith @*/
4577087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
45882bf6240SBarry Smith {
4594ac538c5SBarry Smith   PetscErrorCode ierr;
46082bf6240SBarry Smith 
46182bf6240SBarry Smith   PetscFunctionBegin;
4620700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4634ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
46482bf6240SBarry Smith   PetscFunctionReturn(0);
46582bf6240SBarry Smith }
46682bf6240SBarry Smith 
46782bf6240SBarry Smith 
4684a2ae208SSatish Balay #undef __FUNCT__
4694a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
470ac226902SBarry Smith /*@C
47182bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
47282bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
47382bf6240SBarry Smith 
4743f9fe445SBarry Smith    Logically Collective on TS
47515091d37SBarry Smith 
47682bf6240SBarry Smith    Input Parameters:
47715091d37SBarry Smith +  ts - timestep context
47882bf6240SBarry Smith .  dt - function to compute timestep
47915091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48082bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48182bf6240SBarry Smith 
48215091d37SBarry Smith    Level: intermediate
483fee21e36SBarry Smith 
48482bf6240SBarry Smith    Calling sequence of func:
48587828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
48682bf6240SBarry Smith 
48782bf6240SBarry Smith .  newdt - the newly computed timestep
48882bf6240SBarry Smith .  ctx - [optional] timestep context
48982bf6240SBarry Smith 
49082bf6240SBarry Smith    Notes:
49182bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
49282bf6240SBarry Smith    during the timestepping process.
49382bf6240SBarry Smith 
49482bf6240SBarry Smith .keywords: timestep, pseudo, set
49582bf6240SBarry Smith 
49682bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
49782bf6240SBarry Smith @*/
4987087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
49982bf6240SBarry Smith {
5004ac538c5SBarry Smith   PetscErrorCode ierr;
50182bf6240SBarry Smith 
50282bf6240SBarry Smith   PetscFunctionBegin;
5030700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5044ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
50582bf6240SBarry Smith   PetscFunctionReturn(0);
50682bf6240SBarry Smith }
50782bf6240SBarry Smith 
50882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
50982bf6240SBarry Smith 
510ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
511fb2e594dSBarry Smith EXTERN_C_BEGIN
5124a2ae208SSatish Balay #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5147087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
51582bf6240SBarry Smith {
51682bf6240SBarry Smith   TS_Pseudo *pseudo;
51782bf6240SBarry Smith 
51882bf6240SBarry Smith   PetscFunctionBegin;
51982bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
52082bf6240SBarry Smith   pseudo->verify      = dt;
52182bf6240SBarry Smith   pseudo->verifyctx   = ctx;
52282bf6240SBarry Smith   PetscFunctionReturn(0);
52382bf6240SBarry Smith }
524fb2e594dSBarry Smith EXTERN_C_END
52582bf6240SBarry Smith 
526fb2e594dSBarry Smith EXTERN_C_BEGIN
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5297087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
53082bf6240SBarry Smith {
5314bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53282bf6240SBarry Smith 
53382bf6240SBarry Smith   PetscFunctionBegin;
53482bf6240SBarry Smith   pseudo->dt_increment = inc;
53582bf6240SBarry Smith   PetscFunctionReturn(0);
53682bf6240SBarry Smith }
537fb2e594dSBarry Smith EXTERN_C_END
53882bf6240SBarry Smith 
539fb2e594dSBarry Smith EXTERN_C_BEGIN
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
5427087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
54382bf6240SBarry Smith {
5444bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54582bf6240SBarry Smith 
54682bf6240SBarry Smith   PetscFunctionBegin;
5474bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
54882bf6240SBarry Smith   PetscFunctionReturn(0);
54982bf6240SBarry Smith }
550fb2e594dSBarry Smith EXTERN_C_END
55182bf6240SBarry Smith 
5526849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
553fb2e594dSBarry Smith EXTERN_C_BEGIN
5544a2ae208SSatish Balay #undef __FUNCT__
5554a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
5567087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
55782bf6240SBarry Smith {
5584bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55982bf6240SBarry Smith 
56082bf6240SBarry Smith   PetscFunctionBegin;
56182bf6240SBarry Smith   pseudo->dt      = dt;
56282bf6240SBarry Smith   pseudo->dtctx   = ctx;
56382bf6240SBarry Smith   PetscFunctionReturn(0);
56482bf6240SBarry Smith }
565fb2e594dSBarry Smith EXTERN_C_END
56682bf6240SBarry Smith 
56782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
56810e6a065SJed Brown /*MC
56910e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
57082bf6240SBarry Smith 
57110e6a065SJed Brown   This method solves equations of the form
57210e6a065SJed Brown 
57310e6a065SJed Brown $    F(X,Xdot) = 0
57410e6a065SJed Brown 
57510e6a065SJed Brown   for steady state using the iteration
57610e6a065SJed Brown 
57710e6a065SJed Brown $    [G'] S = -F(X,0)
57810e6a065SJed Brown $    X += S
57910e6a065SJed Brown 
58010e6a065SJed Brown   where
58110e6a065SJed Brown 
58210e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
58310e6a065SJed Brown 
5846f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5856f2d6a7bSJed Brown   state".  See note below.
58610e6a065SJed Brown 
58710e6a065SJed Brown   Options database keys:
58810e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
58910e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
59010e6a065SJed Brown 
59110e6a065SJed Brown   Level: beginner
59210e6a065SJed Brown 
59310e6a065SJed Brown   References:
59410e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
59510e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
59610e6a065SJed Brown 
59710e6a065SJed Brown   Notes:
5986f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
5996f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6006f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6016f2d6a7bSJed Brown 
6026f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6036f2d6a7bSJed Brown 
6046f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6056f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6066f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
60710e6a065SJed Brown 
60810e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
60910e6a065SJed Brown 
61010e6a065SJed Brown M*/
611fb2e594dSBarry Smith EXTERN_C_BEGIN
6124a2ae208SSatish Balay #undef __FUNCT__
6134a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6147087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6152d3f70b5SBarry Smith {
6167bf11e45SBarry Smith   TS_Pseudo      *pseudo;
617dfbe8321SBarry Smith   PetscErrorCode ierr;
6182d3f70b5SBarry Smith 
6193a40ed3dSBarry Smith   PetscFunctionBegin;
620277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
621000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
622000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6232d3f70b5SBarry Smith 
624000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
625000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
626000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6270f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6280f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6297bf11e45SBarry Smith 
6307bf11e45SBarry Smith   /* create the required nonlinear solver context */
631d372ba47SLisandro Dalcin   ts->problem_type = TS_NONLINEAR;
632d372ba47SLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
6332d3f70b5SBarry Smith 
63438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6357bf11e45SBarry Smith   ts->data = (void*)pseudo;
6362d3f70b5SBarry Smith 
63728aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6384bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
63928aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
64082bf6240SBarry Smith 
641f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
642e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6430c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
644f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
645e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6460c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
647f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
648e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6490c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6500c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6510c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6520c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6533a40ed3dSBarry Smith   PetscFunctionReturn(0);
6542d3f70b5SBarry Smith }
655fb2e594dSBarry Smith EXTERN_C_END
6562d3f70b5SBarry Smith 
6574a2ae208SSatish Balay #undef __FUNCT__
6584a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
65982bf6240SBarry Smith /*@C
66082bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
66182bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
66228aa8177SBarry Smith 
66315091d37SBarry Smith    Collective on TS
66415091d37SBarry Smith 
66528aa8177SBarry Smith    Input Parameters:
66628aa8177SBarry Smith .  ts - the timestep context
66782bf6240SBarry Smith .  dtctx - unused timestep context
66828aa8177SBarry Smith 
66982bf6240SBarry Smith    Output Parameter:
67082bf6240SBarry Smith .  newdt - the timestep to use for the next step
67128aa8177SBarry Smith 
67215091d37SBarry Smith    Level: advanced
67315091d37SBarry Smith 
67482bf6240SBarry Smith .keywords: timestep, pseudo, default
675564e8f4eSLois Curfman McInnes 
67682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
67728aa8177SBarry Smith @*/
6787087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
67928aa8177SBarry Smith {
68082bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
68187828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
682dfbe8321SBarry Smith   PetscErrorCode ierr;
68328aa8177SBarry Smith 
6843a40ed3dSBarry Smith   PetscFunctionBegin;
685bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
686*214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
68782bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
68882bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
68982bf6240SBarry Smith     /* first time through so compute initial function norm */
69082bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
69182bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
69282bf6240SBarry Smith   }
69382bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
69482bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
695004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
69682bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
69782bf6240SBarry Smith   } else {
69882bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
69982bf6240SBarry Smith   }
70082bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
70228aa8177SBarry Smith }
703