xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 335f802ec38f19199f77afcc655597b8c823be46)
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;
145c99c2784SJed 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);
207*335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
208*335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr);
209*335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr);
210*335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
211277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
212277b19d0SLisandro Dalcin }
2132d3f70b5SBarry Smith 
2142d3f70b5SBarry Smith /*------------------------------------------------------------*/
2152d3f70b5SBarry Smith 
2164a2ae208SSatish Balay #undef __FUNCT__
2176f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2186f2d6a7bSJed Brown /*
2196f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2206f2d6a7bSJed Brown */
2216f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2222d3f70b5SBarry Smith {
2236f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
224193ac0bcSJed Brown   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
225193ac0bcSJed Brown   PetscScalar    *xdot;
226dfbe8321SBarry Smith   PetscErrorCode ierr;
227a7cc72afSBarry Smith   PetscInt       i,n;
2282d3f70b5SBarry Smith 
2293a40ed3dSBarry Smith   PetscFunctionBegin;
230193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
231193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2326f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2336f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2342d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2356f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2362d3f70b5SBarry Smith   }
237193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
238193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2396f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2406f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2413a40ed3dSBarry Smith   PetscFunctionReturn(0);
2422d3f70b5SBarry Smith }
2432d3f70b5SBarry Smith 
2444a2ae208SSatish Balay #undef __FUNCT__
2450f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2466f2d6a7bSJed Brown /*
2476f2d6a7bSJed Brown     The transient residual is
2486f2d6a7bSJed Brown 
2496f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2506f2d6a7bSJed Brown 
2516f2d6a7bSJed Brown     or for ODE,
2526f2d6a7bSJed Brown 
2536f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2546f2d6a7bSJed Brown 
2556f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2566f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2576f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2586f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2596f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2606f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2616f2d6a7bSJed Brown     residual, and then advances to the next time step.
2626f2d6a7bSJed Brown */
2630f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2642d3f70b5SBarry Smith {
2656f2d6a7bSJed Brown   Vec            Xdot;
266dfbe8321SBarry Smith   PetscErrorCode ierr;
2672d3f70b5SBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionBegin;
2696f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
270193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2716f2d6a7bSJed Brown   PetscFunctionReturn(0);
2726f2d6a7bSJed Brown }
2732d3f70b5SBarry Smith 
2746f2d6a7bSJed Brown #undef __FUNCT__
2750f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2766f2d6a7bSJed Brown /*
2776f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2786f2d6a7bSJed Brown 
2796f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2806f2d6a7bSJed Brown 
2816f2d6a7bSJed Brown     and for ODE:
2826f2d6a7bSJed Brown 
2836f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2846f2d6a7bSJed Brown */
2850f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2866f2d6a7bSJed Brown {
2876f2d6a7bSJed Brown   Vec            Xdot;
2886f2d6a7bSJed Brown   PetscErrorCode ierr;
2896f2d6a7bSJed Brown 
2906f2d6a7bSJed Brown   PetscFunctionBegin;
2916f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
292193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
2933a40ed3dSBarry Smith   PetscFunctionReturn(0);
2942d3f70b5SBarry Smith }
2952d3f70b5SBarry Smith 
2962d3f70b5SBarry Smith 
2974a2ae208SSatish Balay #undef __FUNCT__
2984a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2996849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3002d3f70b5SBarry Smith {
3017bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
302dfbe8321SBarry Smith   PetscErrorCode ierr;
3032d3f70b5SBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
3057bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3067bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3076f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3092d3f70b5SBarry Smith }
3102d3f70b5SBarry Smith /*------------------------------------------------------------*/
3112d3f70b5SBarry Smith 
3124a2ae208SSatish Balay #undef __FUNCT__
313a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
314649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3152d3f70b5SBarry Smith {
3167bf11e45SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
317dfbe8321SBarry Smith   PetscErrorCode   ierr;
318649052a6SBarry Smith   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3192d3f70b5SBarry Smith 
3203a40ed3dSBarry Smith   PetscFunctionBegin;
321193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
322193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
323193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
324193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
325193ac0bcSJed Brown   }
326649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
327649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
328649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3293a40ed3dSBarry Smith   PetscFunctionReturn(0);
3302d3f70b5SBarry Smith }
3312d3f70b5SBarry Smith 
3324a2ae208SSatish Balay #undef __FUNCT__
3334a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3346849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3352d3f70b5SBarry Smith {
3364bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
337dfbe8321SBarry Smith   PetscErrorCode  ierr;
338ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
339649052a6SBarry Smith   PetscViewer     viewer;
3402d3f70b5SBarry Smith 
3413a40ed3dSBarry Smith   PetscFunctionBegin;
342b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
343193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3442d3f70b5SBarry Smith     if (flg) {
345649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
346649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
34728aa8177SBarry Smith     }
34890d69ab7SBarry Smith     flg  = PETSC_FALSE;
349acfcf0e5SJed 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);
350ca4b7087SBarry Smith     if (flg) {
351ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
352ca4b7087SBarry Smith     }
353419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
354b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
3562d3f70b5SBarry Smith }
3572d3f70b5SBarry Smith 
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3606849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3612d3f70b5SBarry Smith {
3623a40ed3dSBarry Smith   PetscFunctionBegin;
3633a40ed3dSBarry Smith   PetscFunctionReturn(0);
3642d3f70b5SBarry Smith }
3652d3f70b5SBarry Smith 
36682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3674a2ae208SSatish Balay #undef __FUNCT__
3684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
369ac226902SBarry Smith /*@C
37082bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
37182bf6240SBarry Smith    last timestep.
37282bf6240SBarry Smith 
3733f9fe445SBarry Smith    Logically Collective on TS
37415091d37SBarry Smith 
37582bf6240SBarry Smith    Input Parameters:
37615091d37SBarry Smith +  ts - timestep context
37782bf6240SBarry Smith .  dt - user-defined function to verify timestep
37815091d37SBarry Smith -  ctx - [optional] user-defined context for private data
37982bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
38082bf6240SBarry Smith 
38115091d37SBarry Smith    Level: advanced
382fee21e36SBarry Smith 
38382bf6240SBarry Smith    Calling sequence of func:
384ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
38582bf6240SBarry Smith 
38682bf6240SBarry Smith .  update - latest solution vector
38782bf6240SBarry Smith .  ctx - [optional] timestep context
38882bf6240SBarry Smith .  newdt - the timestep to use for the next step
38982bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
39082bf6240SBarry Smith 
39182bf6240SBarry Smith    Notes:
39282bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
39382bf6240SBarry Smith    during the timestepping process.
39482bf6240SBarry Smith 
39582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
39682bf6240SBarry Smith 
39782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
39882bf6240SBarry Smith @*/
3997087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
40082bf6240SBarry Smith {
4014ac538c5SBarry Smith   PetscErrorCode ierr;
40282bf6240SBarry Smith 
40382bf6240SBarry Smith   PetscFunctionBegin;
4040700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4054ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
40682bf6240SBarry Smith   PetscFunctionReturn(0);
40782bf6240SBarry Smith }
40882bf6240SBarry Smith 
4094a2ae208SSatish Balay #undef __FUNCT__
4104a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
41182bf6240SBarry Smith /*@
41282bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
41382bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
41482bf6240SBarry Smith 
4153f9fe445SBarry Smith    Logically Collective on TS
416fee21e36SBarry Smith 
41715091d37SBarry Smith     Input Parameters:
41815091d37SBarry Smith +   ts - the timestep context
41915091d37SBarry Smith -   inc - the scaling factor >= 1.0
42015091d37SBarry Smith 
42182bf6240SBarry Smith     Options Database Key:
42282bf6240SBarry Smith $    -ts_pseudo_increment <increment>
42382bf6240SBarry Smith 
42415091d37SBarry Smith     Level: advanced
42515091d37SBarry Smith 
42682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42782bf6240SBarry Smith 
42882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42982bf6240SBarry Smith @*/
4307087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
43182bf6240SBarry Smith {
4324ac538c5SBarry Smith   PetscErrorCode ierr;
43382bf6240SBarry Smith 
43482bf6240SBarry Smith   PetscFunctionBegin;
4350700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
436c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4374ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
43882bf6240SBarry Smith   PetscFunctionReturn(0);
43982bf6240SBarry Smith }
44082bf6240SBarry Smith 
4414a2ae208SSatish Balay #undef __FUNCT__
4424a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
44382bf6240SBarry Smith /*@
44482bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
44582bf6240SBarry Smith     is computed via the formula
44682bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44782bf6240SBarry Smith       rather than the default update,
44882bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44982bf6240SBarry Smith 
4503f9fe445SBarry Smith    Logically Collective on TS
45115091d37SBarry Smith 
45282bf6240SBarry Smith     Input Parameter:
45382bf6240SBarry Smith .   ts - the timestep context
45482bf6240SBarry Smith 
45582bf6240SBarry Smith     Options Database Key:
45682bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45782bf6240SBarry Smith 
45815091d37SBarry Smith     Level: advanced
45915091d37SBarry Smith 
46082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
46182bf6240SBarry Smith 
46282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
46382bf6240SBarry Smith @*/
4647087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
46582bf6240SBarry Smith {
4664ac538c5SBarry Smith   PetscErrorCode ierr;
46782bf6240SBarry Smith 
46882bf6240SBarry Smith   PetscFunctionBegin;
4690700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4704ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
47182bf6240SBarry Smith   PetscFunctionReturn(0);
47282bf6240SBarry Smith }
47382bf6240SBarry Smith 
47482bf6240SBarry Smith 
4754a2ae208SSatish Balay #undef __FUNCT__
4764a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
477ac226902SBarry Smith /*@C
47882bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
47982bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
48082bf6240SBarry Smith 
4813f9fe445SBarry Smith    Logically Collective on TS
48215091d37SBarry Smith 
48382bf6240SBarry Smith    Input Parameters:
48415091d37SBarry Smith +  ts - timestep context
48582bf6240SBarry Smith .  dt - function to compute timestep
48615091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48782bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48882bf6240SBarry Smith 
48915091d37SBarry Smith    Level: intermediate
490fee21e36SBarry Smith 
49182bf6240SBarry Smith    Calling sequence of func:
49287828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
49382bf6240SBarry Smith 
49482bf6240SBarry Smith .  newdt - the newly computed timestep
49582bf6240SBarry Smith .  ctx - [optional] timestep context
49682bf6240SBarry Smith 
49782bf6240SBarry Smith    Notes:
49882bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
49982bf6240SBarry Smith    during the timestepping process.
50082bf6240SBarry Smith 
50182bf6240SBarry Smith .keywords: timestep, pseudo, set
50282bf6240SBarry Smith 
50382bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
50482bf6240SBarry Smith @*/
5057087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
50682bf6240SBarry Smith {
5074ac538c5SBarry Smith   PetscErrorCode ierr;
50882bf6240SBarry Smith 
50982bf6240SBarry Smith   PetscFunctionBegin;
5100700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5114ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
51282bf6240SBarry Smith   PetscFunctionReturn(0);
51382bf6240SBarry Smith }
51482bf6240SBarry Smith 
51582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
51682bf6240SBarry Smith 
517ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
518fb2e594dSBarry Smith EXTERN_C_BEGIN
5194a2ae208SSatish Balay #undef __FUNCT__
5204a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5217087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
52282bf6240SBarry Smith {
52382bf6240SBarry Smith   TS_Pseudo *pseudo;
52482bf6240SBarry Smith 
52582bf6240SBarry Smith   PetscFunctionBegin;
52682bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
52782bf6240SBarry Smith   pseudo->verify      = dt;
52882bf6240SBarry Smith   pseudo->verifyctx   = ctx;
52982bf6240SBarry Smith   PetscFunctionReturn(0);
53082bf6240SBarry Smith }
531fb2e594dSBarry Smith EXTERN_C_END
53282bf6240SBarry Smith 
533fb2e594dSBarry Smith EXTERN_C_BEGIN
5344a2ae208SSatish Balay #undef __FUNCT__
5354a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5367087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
53782bf6240SBarry Smith {
5384bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53982bf6240SBarry Smith 
54082bf6240SBarry Smith   PetscFunctionBegin;
54182bf6240SBarry Smith   pseudo->dt_increment = inc;
54282bf6240SBarry Smith   PetscFunctionReturn(0);
54382bf6240SBarry Smith }
544fb2e594dSBarry Smith EXTERN_C_END
54582bf6240SBarry Smith 
546fb2e594dSBarry Smith EXTERN_C_BEGIN
5474a2ae208SSatish Balay #undef __FUNCT__
5484a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
5497087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
55082bf6240SBarry Smith {
5514bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55282bf6240SBarry Smith 
55382bf6240SBarry Smith   PetscFunctionBegin;
5544bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
55582bf6240SBarry Smith   PetscFunctionReturn(0);
55682bf6240SBarry Smith }
557fb2e594dSBarry Smith EXTERN_C_END
55882bf6240SBarry Smith 
5596849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
560fb2e594dSBarry Smith EXTERN_C_BEGIN
5614a2ae208SSatish Balay #undef __FUNCT__
5624a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
5637087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
56482bf6240SBarry Smith {
5654bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
56682bf6240SBarry Smith 
56782bf6240SBarry Smith   PetscFunctionBegin;
56882bf6240SBarry Smith   pseudo->dt      = dt;
56982bf6240SBarry Smith   pseudo->dtctx   = ctx;
57082bf6240SBarry Smith   PetscFunctionReturn(0);
57182bf6240SBarry Smith }
572fb2e594dSBarry Smith EXTERN_C_END
57382bf6240SBarry Smith 
57482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
57510e6a065SJed Brown /*MC
57610e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
57782bf6240SBarry Smith 
57810e6a065SJed Brown   This method solves equations of the form
57910e6a065SJed Brown 
58010e6a065SJed Brown $    F(X,Xdot) = 0
58110e6a065SJed Brown 
58210e6a065SJed Brown   for steady state using the iteration
58310e6a065SJed Brown 
58410e6a065SJed Brown $    [G'] S = -F(X,0)
58510e6a065SJed Brown $    X += S
58610e6a065SJed Brown 
58710e6a065SJed Brown   where
58810e6a065SJed Brown 
58910e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
59010e6a065SJed Brown 
5916f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5926f2d6a7bSJed Brown   state".  See note below.
59310e6a065SJed Brown 
59410e6a065SJed Brown   Options database keys:
59510e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
59610e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
59710e6a065SJed Brown 
59810e6a065SJed Brown   Level: beginner
59910e6a065SJed Brown 
60010e6a065SJed Brown   References:
60110e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
60210e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
60310e6a065SJed Brown 
60410e6a065SJed Brown   Notes:
6056f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6066f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6076f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6086f2d6a7bSJed Brown 
6096f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6106f2d6a7bSJed Brown 
6116f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6126f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6136f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
61410e6a065SJed Brown 
61510e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
61610e6a065SJed Brown 
61710e6a065SJed Brown M*/
618fb2e594dSBarry Smith EXTERN_C_BEGIN
6194a2ae208SSatish Balay #undef __FUNCT__
6204a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6217087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6222d3f70b5SBarry Smith {
6237bf11e45SBarry Smith   TS_Pseudo      *pseudo;
624dfbe8321SBarry Smith   PetscErrorCode ierr;
625193ac0bcSJed Brown   SNES           snes;
626193ac0bcSJed Brown   const SNESType stype;
6272d3f70b5SBarry Smith 
6283a40ed3dSBarry Smith   PetscFunctionBegin;
629277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
630000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
631000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6322d3f70b5SBarry Smith 
633000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
634000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
635000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6360f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6370f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6387bf11e45SBarry Smith 
639193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
640193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
641193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6422d3f70b5SBarry Smith 
64338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6447bf11e45SBarry Smith   ts->data = (void*)pseudo;
6452d3f70b5SBarry Smith 
64628aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6474bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
64828aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
649193ac0bcSJed Brown   pseudo->fnorm                        = -1;
65082bf6240SBarry Smith 
651f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
652e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6530c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
654f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
655e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6560c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
657f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
658e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6590c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6600c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6610c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6620c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
6642d3f70b5SBarry Smith }
665fb2e594dSBarry Smith EXTERN_C_END
6662d3f70b5SBarry Smith 
6674a2ae208SSatish Balay #undef __FUNCT__
6684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
66982bf6240SBarry Smith /*@C
67082bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
67182bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
67228aa8177SBarry Smith 
67315091d37SBarry Smith    Collective on TS
67415091d37SBarry Smith 
67528aa8177SBarry Smith    Input Parameters:
67628aa8177SBarry Smith .  ts - the timestep context
67782bf6240SBarry Smith .  dtctx - unused timestep context
67828aa8177SBarry Smith 
67982bf6240SBarry Smith    Output Parameter:
68082bf6240SBarry Smith .  newdt - the timestep to use for the next step
68128aa8177SBarry Smith 
68215091d37SBarry Smith    Level: advanced
68315091d37SBarry Smith 
68482bf6240SBarry Smith .keywords: timestep, pseudo, default
685564e8f4eSLois Curfman McInnes 
68682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
68728aa8177SBarry Smith @*/
6887087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
68928aa8177SBarry Smith {
69082bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
69187828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
692dfbe8321SBarry Smith   PetscErrorCode ierr;
69328aa8177SBarry Smith 
6943a40ed3dSBarry Smith   PetscFunctionBegin;
695bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
696214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
69782bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
69882bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
69982bf6240SBarry Smith     /* first time through so compute initial function norm */
70082bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
70182bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
70282bf6240SBarry Smith   }
70382bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
70482bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
705004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
70682bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
70782bf6240SBarry Smith   } else {
70882bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
70982bf6240SBarry Smith   }
71082bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7113a40ed3dSBarry Smith   PetscFunctionReturn(0);
71228aa8177SBarry Smith }
713