xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 193ac0bc5128976c4ec2c1a67f3f9cb026b77f22)
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"
141*193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1422d3f70b5SBarry Smith {
143277b19d0SLisandro Dalcin   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
144dfbe8321SBarry Smith   PetscErrorCode ierr;
145*193ac0bcSJed Brown   PetscInt       i,its,lits,rej;
146ace3abfcSBarry Smith   PetscBool      ok;
147*193ac0bcSJed Brown   PetscReal      new_time_step;
148*193ac0bcSJed Brown   SNESConvergedReason snesreason;
1492d3f70b5SBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1513f2090d5SJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
152*193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
153*193ac0bcSJed Brown   ierr = TSPseudoComputeTimeStep(ts,&ts->next_time_step);CHKERRQ(ierr);
154*193ac0bcSJed Brown   for (rej=0; ; ) {
155*193ac0bcSJed Brown     ts->time_step = ts->next_time_step;
156*193ac0bcSJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr);
157b850b91aSLisandro Dalcin     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1589f954729SBarry Smith     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
159c9780b6fSBarry Smith     ts->nonlinear_its += its; ts->linear_its += lits;
160*193ac0bcSJed Brown     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
161*193ac0bcSJed Brown     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
162*193ac0bcSJed Brown     if (snesreason < 0) {
163*193ac0bcSJed Brown       if (++ts->num_snes_failures >= ts->max_snes_failures) {
164*193ac0bcSJed Brown         ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
165*193ac0bcSJed 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);
166*193ac0bcSJed Brown         break;
1677bf11e45SBarry Smith       }
168*193ac0bcSJed Brown     }
169*193ac0bcSJed Brown     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->next_time_step,&ok);CHKERRQ(ierr);
170*193ac0bcSJed Brown     printf("step %d, ok %d\n",rej,ok);
171*193ac0bcSJed Brown     if (ok) {
172*193ac0bcSJed Brown       ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
173*193ac0bcSJed Brown       break;
174*193ac0bcSJed Brown     }
175*193ac0bcSJed Brown     ts->reject++;
176*193ac0bcSJed Brown     if (++rej >= ts->max_reject) {
177*193ac0bcSJed Brown       ts->reason = TS_DIVERGED_STEP_REJECTED;
178*193ac0bcSJed Brown     }
179*193ac0bcSJed Brown   }
180*193ac0bcSJed Brown   ts->ptime += ts->time_step;
1812d3f70b5SBarry Smith   ts->steps++;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1832d3f70b5SBarry Smith }
1842d3f70b5SBarry Smith 
1852d3f70b5SBarry Smith /*------------------------------------------------------------*/
1864a2ae208SSatish Balay #undef __FUNCT__
187277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
188277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1892d3f70b5SBarry Smith {
1907bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
191dfbe8321SBarry Smith   PetscErrorCode ierr;
1922d3f70b5SBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
1946bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
1956bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
1966bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1982d3f70b5SBarry Smith }
1992d3f70b5SBarry Smith 
200277b19d0SLisandro Dalcin #undef __FUNCT__
201277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
202277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
203277b19d0SLisandro Dalcin {
204277b19d0SLisandro Dalcin   PetscErrorCode ierr;
205277b19d0SLisandro Dalcin 
206277b19d0SLisandro Dalcin   PetscFunctionBegin;
207277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
208277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
209277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
210277b19d0SLisandro Dalcin }
2112d3f70b5SBarry Smith 
2122d3f70b5SBarry Smith /*------------------------------------------------------------*/
2132d3f70b5SBarry Smith 
2144a2ae208SSatish Balay #undef __FUNCT__
2156f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2166f2d6a7bSJed Brown /*
2176f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2186f2d6a7bSJed Brown */
2196f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2202d3f70b5SBarry Smith {
2216f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
222*193ac0bcSJed Brown   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
223*193ac0bcSJed Brown   PetscScalar    *xdot;
224dfbe8321SBarry Smith   PetscErrorCode ierr;
225a7cc72afSBarry Smith   PetscInt       i,n;
2262d3f70b5SBarry Smith 
2273a40ed3dSBarry Smith   PetscFunctionBegin;
228*193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
229*193ac0bcSJed Brown   ierr = VecGetArrayRead(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   }
235*193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
236*193ac0bcSJed Brown   ierr = VecRestoreArrayRead(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*193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,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*193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,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;
319*193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
320*193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
321*193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
322*193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
323*193ac0bcSJed Brown   }
324649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
325649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
326649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3273a40ed3dSBarry Smith   PetscFunctionReturn(0);
3282d3f70b5SBarry Smith }
3292d3f70b5SBarry Smith 
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3326849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3332d3f70b5SBarry Smith {
3344bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
335dfbe8321SBarry Smith   PetscErrorCode  ierr;
336ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
337649052a6SBarry Smith   PetscViewer     viewer;
3382d3f70b5SBarry Smith 
3393a40ed3dSBarry Smith   PetscFunctionBegin;
340b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
341*193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3422d3f70b5SBarry Smith     if (flg) {
343649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
344649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
34528aa8177SBarry Smith     }
34690d69ab7SBarry Smith     flg  = PETSC_FALSE;
347acfcf0e5SJed 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);
348ca4b7087SBarry Smith     if (flg) {
349ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
350ca4b7087SBarry Smith     }
351419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
352b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
3542d3f70b5SBarry Smith }
3552d3f70b5SBarry Smith 
3564a2ae208SSatish Balay #undef __FUNCT__
3574a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3586849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3592d3f70b5SBarry Smith {
3603a40ed3dSBarry Smith   PetscFunctionBegin;
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
3622d3f70b5SBarry Smith }
3632d3f70b5SBarry Smith 
36482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3654a2ae208SSatish Balay #undef __FUNCT__
3664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
367ac226902SBarry Smith /*@C
36882bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
36982bf6240SBarry Smith    last timestep.
37082bf6240SBarry Smith 
3713f9fe445SBarry Smith    Logically Collective on TS
37215091d37SBarry Smith 
37382bf6240SBarry Smith    Input Parameters:
37415091d37SBarry Smith +  ts - timestep context
37582bf6240SBarry Smith .  dt - user-defined function to verify timestep
37615091d37SBarry Smith -  ctx - [optional] user-defined context for private data
37782bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
37882bf6240SBarry Smith 
37915091d37SBarry Smith    Level: advanced
380fee21e36SBarry Smith 
38182bf6240SBarry Smith    Calling sequence of func:
382ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
38382bf6240SBarry Smith 
38482bf6240SBarry Smith .  update - latest solution vector
38582bf6240SBarry Smith .  ctx - [optional] timestep context
38682bf6240SBarry Smith .  newdt - the timestep to use for the next step
38782bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
38882bf6240SBarry Smith 
38982bf6240SBarry Smith    Notes:
39082bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
39182bf6240SBarry Smith    during the timestepping process.
39282bf6240SBarry Smith 
39382bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
39482bf6240SBarry Smith 
39582bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
39682bf6240SBarry Smith @*/
3977087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
39882bf6240SBarry Smith {
3994ac538c5SBarry Smith   PetscErrorCode ierr;
40082bf6240SBarry Smith 
40182bf6240SBarry Smith   PetscFunctionBegin;
4020700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4034ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
40482bf6240SBarry Smith   PetscFunctionReturn(0);
40582bf6240SBarry Smith }
40682bf6240SBarry Smith 
4074a2ae208SSatish Balay #undef __FUNCT__
4084a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
40982bf6240SBarry Smith /*@
41082bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
41182bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
41282bf6240SBarry Smith 
4133f9fe445SBarry Smith    Logically Collective on TS
414fee21e36SBarry Smith 
41515091d37SBarry Smith     Input Parameters:
41615091d37SBarry Smith +   ts - the timestep context
41715091d37SBarry Smith -   inc - the scaling factor >= 1.0
41815091d37SBarry Smith 
41982bf6240SBarry Smith     Options Database Key:
42082bf6240SBarry Smith $    -ts_pseudo_increment <increment>
42182bf6240SBarry Smith 
42215091d37SBarry Smith     Level: advanced
42315091d37SBarry Smith 
42482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42582bf6240SBarry Smith 
42682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42782bf6240SBarry Smith @*/
4287087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
42982bf6240SBarry Smith {
4304ac538c5SBarry Smith   PetscErrorCode ierr;
43182bf6240SBarry Smith 
43282bf6240SBarry Smith   PetscFunctionBegin;
4330700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
434c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4354ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
43682bf6240SBarry Smith   PetscFunctionReturn(0);
43782bf6240SBarry Smith }
43882bf6240SBarry Smith 
4394a2ae208SSatish Balay #undef __FUNCT__
4404a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
44182bf6240SBarry Smith /*@
44282bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
44382bf6240SBarry Smith     is computed via the formula
44482bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44582bf6240SBarry Smith       rather than the default update,
44682bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44782bf6240SBarry Smith 
4483f9fe445SBarry Smith    Logically Collective on TS
44915091d37SBarry Smith 
45082bf6240SBarry Smith     Input Parameter:
45182bf6240SBarry Smith .   ts - the timestep context
45282bf6240SBarry Smith 
45382bf6240SBarry Smith     Options Database Key:
45482bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45582bf6240SBarry Smith 
45615091d37SBarry Smith     Level: advanced
45715091d37SBarry Smith 
45882bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
45982bf6240SBarry Smith 
46082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
46182bf6240SBarry Smith @*/
4627087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
46382bf6240SBarry Smith {
4644ac538c5SBarry Smith   PetscErrorCode ierr;
46582bf6240SBarry Smith 
46682bf6240SBarry Smith   PetscFunctionBegin;
4670700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4684ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
46982bf6240SBarry Smith   PetscFunctionReturn(0);
47082bf6240SBarry Smith }
47182bf6240SBarry Smith 
47282bf6240SBarry Smith 
4734a2ae208SSatish Balay #undef __FUNCT__
4744a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
475ac226902SBarry Smith /*@C
47682bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
47782bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
47882bf6240SBarry Smith 
4793f9fe445SBarry Smith    Logically Collective on TS
48015091d37SBarry Smith 
48182bf6240SBarry Smith    Input Parameters:
48215091d37SBarry Smith +  ts - timestep context
48382bf6240SBarry Smith .  dt - function to compute timestep
48415091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48582bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48682bf6240SBarry Smith 
48715091d37SBarry Smith    Level: intermediate
488fee21e36SBarry Smith 
48982bf6240SBarry Smith    Calling sequence of func:
49087828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
49182bf6240SBarry Smith 
49282bf6240SBarry Smith .  newdt - the newly computed timestep
49382bf6240SBarry Smith .  ctx - [optional] timestep context
49482bf6240SBarry Smith 
49582bf6240SBarry Smith    Notes:
49682bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
49782bf6240SBarry Smith    during the timestepping process.
49882bf6240SBarry Smith 
49982bf6240SBarry Smith .keywords: timestep, pseudo, set
50082bf6240SBarry Smith 
50182bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
50282bf6240SBarry Smith @*/
5037087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
50482bf6240SBarry Smith {
5054ac538c5SBarry Smith   PetscErrorCode ierr;
50682bf6240SBarry Smith 
50782bf6240SBarry Smith   PetscFunctionBegin;
5080700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5094ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
51082bf6240SBarry Smith   PetscFunctionReturn(0);
51182bf6240SBarry Smith }
51282bf6240SBarry Smith 
51382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
51482bf6240SBarry Smith 
515ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
516fb2e594dSBarry Smith EXTERN_C_BEGIN
5174a2ae208SSatish Balay #undef __FUNCT__
5184a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5197087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
52082bf6240SBarry Smith {
52182bf6240SBarry Smith   TS_Pseudo *pseudo;
52282bf6240SBarry Smith 
52382bf6240SBarry Smith   PetscFunctionBegin;
52482bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
52582bf6240SBarry Smith   pseudo->verify      = dt;
52682bf6240SBarry Smith   pseudo->verifyctx   = ctx;
52782bf6240SBarry Smith   PetscFunctionReturn(0);
52882bf6240SBarry Smith }
529fb2e594dSBarry Smith EXTERN_C_END
53082bf6240SBarry Smith 
531fb2e594dSBarry Smith EXTERN_C_BEGIN
5324a2ae208SSatish Balay #undef __FUNCT__
5334a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5347087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
53582bf6240SBarry Smith {
5364bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53782bf6240SBarry Smith 
53882bf6240SBarry Smith   PetscFunctionBegin;
53982bf6240SBarry Smith   pseudo->dt_increment = inc;
54082bf6240SBarry Smith   PetscFunctionReturn(0);
54182bf6240SBarry Smith }
542fb2e594dSBarry Smith EXTERN_C_END
54382bf6240SBarry Smith 
544fb2e594dSBarry Smith EXTERN_C_BEGIN
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
5477087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
54882bf6240SBarry Smith {
5494bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55082bf6240SBarry Smith 
55182bf6240SBarry Smith   PetscFunctionBegin;
5524bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
55382bf6240SBarry Smith   PetscFunctionReturn(0);
55482bf6240SBarry Smith }
555fb2e594dSBarry Smith EXTERN_C_END
55682bf6240SBarry Smith 
5576849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
558fb2e594dSBarry Smith EXTERN_C_BEGIN
5594a2ae208SSatish Balay #undef __FUNCT__
5604a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
5617087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
56282bf6240SBarry Smith {
5634bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
56482bf6240SBarry Smith 
56582bf6240SBarry Smith   PetscFunctionBegin;
56682bf6240SBarry Smith   pseudo->dt      = dt;
56782bf6240SBarry Smith   pseudo->dtctx   = ctx;
56882bf6240SBarry Smith   PetscFunctionReturn(0);
56982bf6240SBarry Smith }
570fb2e594dSBarry Smith EXTERN_C_END
57182bf6240SBarry Smith 
57282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
57310e6a065SJed Brown /*MC
57410e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
57582bf6240SBarry Smith 
57610e6a065SJed Brown   This method solves equations of the form
57710e6a065SJed Brown 
57810e6a065SJed Brown $    F(X,Xdot) = 0
57910e6a065SJed Brown 
58010e6a065SJed Brown   for steady state using the iteration
58110e6a065SJed Brown 
58210e6a065SJed Brown $    [G'] S = -F(X,0)
58310e6a065SJed Brown $    X += S
58410e6a065SJed Brown 
58510e6a065SJed Brown   where
58610e6a065SJed Brown 
58710e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
58810e6a065SJed Brown 
5896f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5906f2d6a7bSJed Brown   state".  See note below.
59110e6a065SJed Brown 
59210e6a065SJed Brown   Options database keys:
59310e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
59410e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
59510e6a065SJed Brown 
59610e6a065SJed Brown   Level: beginner
59710e6a065SJed Brown 
59810e6a065SJed Brown   References:
59910e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
60010e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
60110e6a065SJed Brown 
60210e6a065SJed Brown   Notes:
6036f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6046f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6056f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6066f2d6a7bSJed Brown 
6076f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6086f2d6a7bSJed Brown 
6096f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6106f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6116f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
61210e6a065SJed Brown 
61310e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
61410e6a065SJed Brown 
61510e6a065SJed Brown M*/
616fb2e594dSBarry Smith EXTERN_C_BEGIN
6174a2ae208SSatish Balay #undef __FUNCT__
6184a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6197087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6202d3f70b5SBarry Smith {
6217bf11e45SBarry Smith   TS_Pseudo      *pseudo;
622dfbe8321SBarry Smith   PetscErrorCode ierr;
623*193ac0bcSJed Brown   SNES           snes;
624*193ac0bcSJed Brown   const SNESType stype;
6252d3f70b5SBarry Smith 
6263a40ed3dSBarry Smith   PetscFunctionBegin;
627277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
628000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
629000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6302d3f70b5SBarry Smith 
631000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
632000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
633000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6340f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6350f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6367bf11e45SBarry Smith 
637*193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
638*193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
639*193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6402d3f70b5SBarry Smith 
64138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6427bf11e45SBarry Smith   ts->data = (void*)pseudo;
6432d3f70b5SBarry Smith 
64428aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6454bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
64628aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
647*193ac0bcSJed Brown   pseudo->fnorm                        = -1;
64882bf6240SBarry Smith 
649f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
650e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6510c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
652f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
653e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6540c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
655f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
656e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6570c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6580c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6590c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6600c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
6622d3f70b5SBarry Smith }
663fb2e594dSBarry Smith EXTERN_C_END
6642d3f70b5SBarry Smith 
6654a2ae208SSatish Balay #undef __FUNCT__
6664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
66782bf6240SBarry Smith /*@C
66882bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
66982bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
67028aa8177SBarry Smith 
67115091d37SBarry Smith    Collective on TS
67215091d37SBarry Smith 
67328aa8177SBarry Smith    Input Parameters:
67428aa8177SBarry Smith .  ts - the timestep context
67582bf6240SBarry Smith .  dtctx - unused timestep context
67628aa8177SBarry Smith 
67782bf6240SBarry Smith    Output Parameter:
67882bf6240SBarry Smith .  newdt - the timestep to use for the next step
67928aa8177SBarry Smith 
68015091d37SBarry Smith    Level: advanced
68115091d37SBarry Smith 
68282bf6240SBarry Smith .keywords: timestep, pseudo, default
683564e8f4eSLois Curfman McInnes 
68482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
68528aa8177SBarry Smith @*/
6867087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
68728aa8177SBarry Smith {
68882bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
68987828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
690dfbe8321SBarry Smith   PetscErrorCode ierr;
69128aa8177SBarry Smith 
6923a40ed3dSBarry Smith   PetscFunctionBegin;
693bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
694214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
69582bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
69682bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
69782bf6240SBarry Smith     /* first time through so compute initial function norm */
69882bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
69982bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
70082bf6240SBarry Smith   }
70182bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
70282bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
703004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
70482bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
70582bf6240SBarry Smith   } else {
70682bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
70782bf6240SBarry Smith   }
70882bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7093a40ed3dSBarry Smith   PetscFunctionReturn(0);
71028aa8177SBarry Smith }
711