xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
32d3f70b5SBarry Smith /*
4fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
52d3f70b5SBarry Smith */
67c4f633dSBarry Smith #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
72d3f70b5SBarry Smith 
82d3f70b5SBarry Smith typedef struct {
92d3f70b5SBarry Smith   Vec  update;      /* work vector where new solution is formed */
102d3f70b5SBarry Smith   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
116f2d6a7bSJed Brown   Vec  xdot;        /* work vector for time derivative of state */
122d3f70b5SBarry Smith 
132d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
142d3f70b5SBarry Smith 
156849ba73SBarry Smith   PetscErrorCode (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
162d3f70b5SBarry Smith   void           *dtctx;
17a7cc72afSBarry Smith   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscTruth*); /* verify previous timestep and related context */
187bf11e45SBarry Smith   void           *verifyctx;
192d3f70b5SBarry Smith 
2087828ca2SBarry Smith   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
2187828ca2SBarry Smith   PetscReal  fnorm_previous;
2228aa8177SBarry Smith 
2387828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
244bbc92c1SBarry Smith   PetscTruth increment_dt_from_initial_dt;
257bf11e45SBarry Smith } TS_Pseudo;
262d3f70b5SBarry Smith 
272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
282d3f70b5SBarry Smith 
294a2ae208SSatish Balay #undef __FUNCT__
304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
317bf11e45SBarry Smith /*@
327bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
33564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
342d3f70b5SBarry Smith 
3515091d37SBarry Smith     Collective on TS
3615091d37SBarry Smith 
377bf11e45SBarry Smith     Input Parameter:
387bf11e45SBarry Smith .   ts - timestep context
397bf11e45SBarry Smith 
407bf11e45SBarry Smith     Output Parameter:
41fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
42fb4a63b6SLois Curfman McInnes 
4315091d37SBarry Smith     Level: advanced
44564e8f4eSLois Curfman McInnes 
45564e8f4eSLois Curfman McInnes     Notes:
46564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
47564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
48564e8f4eSLois Curfman McInnes 
49fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
50564e8f4eSLois Curfman McInnes 
51564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
527bf11e45SBarry Smith @*/
5363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
547bf11e45SBarry Smith {
557bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
56dfbe8321SBarry Smith   PetscErrorCode ierr;
577bf11e45SBarry Smith 
583a40ed3dSBarry Smith   PetscFunctionBegin;
59d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
607bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
61d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
623a40ed3dSBarry Smith   PetscFunctionReturn(0);
637bf11e45SBarry Smith }
647bf11e45SBarry Smith 
657bf11e45SBarry Smith 
667bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
697bf11e45SBarry Smith /*@C
70639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
717bf11e45SBarry Smith 
7215091d37SBarry Smith    Collective on TS
7315091d37SBarry Smith 
747bf11e45SBarry Smith    Input Parameters:
7515091d37SBarry Smith +  ts - the timestep context
767bf11e45SBarry Smith .  dtctx - unused timestep context
7715091d37SBarry Smith -  update - latest solution vector
787bf11e45SBarry Smith 
79564e8f4eSLois Curfman McInnes    Output Parameters:
8015091d37SBarry Smith +  newdt - the timestep to use for the next step
8115091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
827bf11e45SBarry Smith 
8315091d37SBarry Smith    Level: advanced
84fee21e36SBarry Smith 
85564e8f4eSLois Curfman McInnes    Note:
86564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
87564e8f4eSLois Curfman McInnes    timestep.
88564e8f4eSLois Curfman McInnes 
89564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
90564e8f4eSLois Curfman McInnes 
91564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
927bf11e45SBarry Smith @*/
9363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *flag)
947bf11e45SBarry Smith {
953a40ed3dSBarry Smith   PetscFunctionBegin;
96a7cc72afSBarry Smith   *flag = PETSC_TRUE;
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
987bf11e45SBarry Smith }
997bf11e45SBarry Smith 
1007bf11e45SBarry Smith 
1014a2ae208SSatish Balay #undef __FUNCT__
1024a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1037bf11e45SBarry Smith /*@
104564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1057bf11e45SBarry Smith 
10615091d37SBarry Smith     Collective on TS
10715091d37SBarry Smith 
108fb4a63b6SLois Curfman McInnes     Input Parameters:
10915091d37SBarry Smith +   ts - timestep context
11015091d37SBarry Smith -   update - latest solution vector
1117bf11e45SBarry Smith 
112fb4a63b6SLois Curfman McInnes     Output Parameters:
11315091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11415091d37SBarry Smith -   flag - indicates if current timestep was ok
1157bf11e45SBarry Smith 
11615091d37SBarry Smith     Level: advanced
117fee21e36SBarry Smith 
118564e8f4eSLois Curfman McInnes     Notes:
119564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
120564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
121564e8f4eSLois Curfman McInnes 
122564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
123564e8f4eSLois Curfman McInnes 
124564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1257bf11e45SBarry Smith @*/
12663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *flag)
1277bf11e45SBarry Smith {
1287bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
129dfbe8321SBarry Smith   PetscErrorCode ierr;
1307bf11e45SBarry Smith 
1313a40ed3dSBarry Smith   PetscFunctionBegin;
132a7cc72afSBarry Smith   if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);}
1337bf11e45SBarry Smith 
1347bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1357bf11e45SBarry Smith 
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
1377bf11e45SBarry Smith }
1387bf11e45SBarry Smith 
1397bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1407bf11e45SBarry Smith 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
143a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime)
1442d3f70b5SBarry Smith {
1452d3f70b5SBarry Smith   Vec            sol = ts->vec_sol;
146dfbe8321SBarry Smith   PetscErrorCode ierr;
147a7cc72afSBarry Smith   PetscInt       i,max_steps = ts->max_steps,its,lits;
148a7cc72afSBarry Smith   PetscTruth     ok;
1497bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
15087828ca2SBarry Smith   PetscReal      current_time_step;
1512d3f70b5SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
1532d3f70b5SBarry Smith   *steps = -ts->steps;
1542d3f70b5SBarry Smith 
1557bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
1562d3f70b5SBarry Smith   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
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;
163f69a0ea3SMatthew Knepley       ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->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;
1677bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->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     }
1727bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
1732d3f70b5SBarry Smith     ts->steps++;
1743f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
1752d3f70b5SBarry Smith   }
176e6e5fe25SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
177e6e5fe25SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
178e6e5fe25SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1792d3f70b5SBarry Smith 
1802d3f70b5SBarry Smith   *steps += ts->steps;
181142b95e3SSatish Balay   *ptime  = ts->ptime;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1832d3f70b5SBarry Smith }
1842d3f70b5SBarry Smith 
1852d3f70b5SBarry Smith /*------------------------------------------------------------*/
1864a2ae208SSatish Balay #undef __FUNCT__
1874a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo"
1886849ba73SBarry Smith static PetscErrorCode TSDestroy_Pseudo(TS ts)
1892d3f70b5SBarry Smith {
1907bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
191dfbe8321SBarry Smith   PetscErrorCode ierr;
1922d3f70b5SBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
194e4ef1970SSatish Balay   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
1957bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
1966f2d6a7bSJed Brown   if (pseudo->xdot) {ierr = VecDestroy(pseudo->xdot);CHKERRQ(ierr);}
197606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1992d3f70b5SBarry Smith }
2002d3f70b5SBarry Smith 
2012d3f70b5SBarry Smith 
2022d3f70b5SBarry Smith /*------------------------------------------------------------*/
2032d3f70b5SBarry Smith 
2044a2ae208SSatish Balay #undef __FUNCT__
2056f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2066f2d6a7bSJed Brown /*
2076f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2086f2d6a7bSJed Brown */
2096f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2102d3f70b5SBarry Smith {
2116f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
2126f2d6a7bSJed Brown   PetscScalar    mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot;
213dfbe8321SBarry Smith   PetscErrorCode ierr;
214a7cc72afSBarry Smith   PetscInt       i,n;
2152d3f70b5SBarry Smith 
2163a40ed3dSBarry Smith   PetscFunctionBegin;
2176f2d6a7bSJed Brown   ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2186f2d6a7bSJed Brown   ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr);
2196f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2206f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2212d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2226f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2232d3f70b5SBarry Smith   }
2246f2d6a7bSJed Brown   ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2256f2d6a7bSJed Brown   ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr);
2266f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2276f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
2292d3f70b5SBarry Smith }
2302d3f70b5SBarry Smith 
2314a2ae208SSatish Balay #undef __FUNCT__
2326f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoFunction"
2336f2d6a7bSJed Brown /*
2346f2d6a7bSJed Brown     The transient residual is
2356f2d6a7bSJed Brown 
2366f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2376f2d6a7bSJed Brown 
2386f2d6a7bSJed Brown     or for ODE,
2396f2d6a7bSJed Brown 
2406f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2416f2d6a7bSJed Brown 
2426f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2436f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2446f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2456f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2466f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2476f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2486f2d6a7bSJed Brown     residual, and then advances to the next time step.
2496f2d6a7bSJed Brown */
2506f2d6a7bSJed Brown static PetscErrorCode TSPseudoFunction(SNES snes,Vec X,Vec Y,void *ctx)
2512d3f70b5SBarry Smith {
2522d3f70b5SBarry Smith   TS             ts = (TS)ctx;
2536f2d6a7bSJed Brown   Vec            Xdot;
254dfbe8321SBarry Smith   PetscErrorCode ierr;
2552d3f70b5SBarry Smith 
2563a40ed3dSBarry Smith   PetscFunctionBegin;
2576f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
2586f2d6a7bSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr);
2596f2d6a7bSJed Brown   PetscFunctionReturn(0);
2606f2d6a7bSJed Brown }
2612d3f70b5SBarry Smith 
2626f2d6a7bSJed Brown #undef __FUNCT__
2636f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoJacobian"
2646f2d6a7bSJed Brown /*
2656f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2666f2d6a7bSJed Brown 
2676f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2686f2d6a7bSJed Brown 
2696f2d6a7bSJed Brown     and for ODE:
2706f2d6a7bSJed Brown 
2716f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2726f2d6a7bSJed Brown */
2736f2d6a7bSJed Brown static PetscErrorCode TSPseudoJacobian(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2746f2d6a7bSJed Brown {
2756f2d6a7bSJed Brown   TS             ts = (TS)ctx;
2766f2d6a7bSJed Brown   Vec            Xdot;
2776f2d6a7bSJed Brown   PetscErrorCode ierr;
2786f2d6a7bSJed Brown 
2796f2d6a7bSJed Brown   PetscFunctionBegin;
2806f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
2816f2d6a7bSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr);
2823a40ed3dSBarry Smith   PetscFunctionReturn(0);
2832d3f70b5SBarry Smith }
2842d3f70b5SBarry Smith 
2852d3f70b5SBarry Smith 
2864a2ae208SSatish Balay #undef __FUNCT__
2874a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2886849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2892d3f70b5SBarry Smith {
2907bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
291dfbe8321SBarry Smith   PetscErrorCode ierr;
2922d3f70b5SBarry Smith 
2933a40ed3dSBarry Smith   PetscFunctionBegin;
2947bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2957bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2966f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
2977bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2988beb423aSHong Zhang   ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
3002d3f70b5SBarry Smith }
3012d3f70b5SBarry Smith /*------------------------------------------------------------*/
3022d3f70b5SBarry Smith 
3034a2ae208SSatish Balay #undef __FUNCT__
304a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
305a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
3062d3f70b5SBarry Smith {
3077bf11e45SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
308dfbe8321SBarry Smith   PetscErrorCode          ierr;
309a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
3102d3f70b5SBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionBegin;
312a34d58ebSBarry Smith   if (!ctx) {
3137adad957SLisandro Dalcin     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
314a34d58ebSBarry Smith   }
315a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
316a34d58ebSBarry Smith   if (!ctx) {
317a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
318a34d58ebSBarry Smith   }
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
3202d3f70b5SBarry Smith }
3212d3f70b5SBarry Smith 
3224a2ae208SSatish Balay #undef __FUNCT__
3234a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3246849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3252d3f70b5SBarry Smith {
3264bbc92c1SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
327dfbe8321SBarry Smith   PetscErrorCode          ierr;
32890d69ab7SBarry Smith   PetscTruth              flg = PETSC_FALSE;
329a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer;
3302d3f70b5SBarry Smith 
3313a40ed3dSBarry Smith   PetscFunctionBegin;
332b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
33390d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3342d3f70b5SBarry Smith     if (flg) {
3357adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
336a34d58ebSBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
33728aa8177SBarry Smith     }
33890d69ab7SBarry Smith     flg  = PETSC_FALSE;
33990d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
340ca4b7087SBarry Smith     if (flg) {
341ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
342ca4b7087SBarry Smith     }
343419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
344b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3453a40ed3dSBarry Smith   PetscFunctionReturn(0);
3462d3f70b5SBarry Smith }
3472d3f70b5SBarry Smith 
3484a2ae208SSatish Balay #undef __FUNCT__
3494a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3506849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3512d3f70b5SBarry Smith {
3523a40ed3dSBarry Smith   PetscFunctionBegin;
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
3542d3f70b5SBarry Smith }
3552d3f70b5SBarry Smith 
35682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3574a2ae208SSatish Balay #undef __FUNCT__
3584a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
359ac226902SBarry Smith /*@C
36082bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
36182bf6240SBarry Smith    last timestep.
36282bf6240SBarry Smith 
36315091d37SBarry Smith    Collective on TS
36415091d37SBarry Smith 
36582bf6240SBarry Smith    Input Parameters:
36615091d37SBarry Smith +  ts - timestep context
36782bf6240SBarry Smith .  dt - user-defined function to verify timestep
36815091d37SBarry Smith -  ctx - [optional] user-defined context for private data
36982bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
37082bf6240SBarry Smith 
37115091d37SBarry Smith    Level: advanced
372fee21e36SBarry Smith 
37382bf6240SBarry Smith    Calling sequence of func:
374a7cc72afSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag);
37582bf6240SBarry Smith 
37682bf6240SBarry Smith .  update - latest solution vector
37782bf6240SBarry Smith .  ctx - [optional] timestep context
37882bf6240SBarry Smith .  newdt - the timestep to use for the next step
37982bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
38082bf6240SBarry Smith 
38182bf6240SBarry Smith    Notes:
38282bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
38382bf6240SBarry Smith    during the timestepping process.
38482bf6240SBarry Smith 
38582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
38682bf6240SBarry Smith 
38782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
38882bf6240SBarry Smith @*/
38963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx)
39082bf6240SBarry Smith {
391a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *);
39282bf6240SBarry Smith 
39382bf6240SBarry Smith   PetscFunctionBegin;
394*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
39582bf6240SBarry Smith 
396c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
39782bf6240SBarry Smith   if (f) {
39882bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
39982bf6240SBarry Smith   }
40082bf6240SBarry Smith   PetscFunctionReturn(0);
40182bf6240SBarry Smith }
40282bf6240SBarry Smith 
4034a2ae208SSatish Balay #undef __FUNCT__
4044a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
40582bf6240SBarry Smith /*@
40682bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
40782bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
40882bf6240SBarry Smith 
409fee21e36SBarry Smith    Collective on TS
410fee21e36SBarry Smith 
41115091d37SBarry Smith     Input Parameters:
41215091d37SBarry Smith +   ts - the timestep context
41315091d37SBarry Smith -   inc - the scaling factor >= 1.0
41415091d37SBarry Smith 
41582bf6240SBarry Smith     Options Database Key:
41682bf6240SBarry Smith $    -ts_pseudo_increment <increment>
41782bf6240SBarry Smith 
41815091d37SBarry Smith     Level: advanced
41915091d37SBarry Smith 
42082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42182bf6240SBarry Smith 
42282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42382bf6240SBarry Smith @*/
42463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
42582bf6240SBarry Smith {
426dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscReal);
42782bf6240SBarry Smith 
42882bf6240SBarry Smith   PetscFunctionBegin;
429*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
43082bf6240SBarry Smith 
431c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
43282bf6240SBarry Smith   if (f) {
43382bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
43482bf6240SBarry Smith   }
43582bf6240SBarry Smith   PetscFunctionReturn(0);
43682bf6240SBarry Smith }
43782bf6240SBarry Smith 
4384a2ae208SSatish Balay #undef __FUNCT__
4394a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
44082bf6240SBarry Smith /*@
44182bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
44282bf6240SBarry Smith     is computed via the formula
44382bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44482bf6240SBarry Smith       rather than the default update,
44582bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44682bf6240SBarry Smith 
44715091d37SBarry Smith    Collective on TS
44815091d37SBarry Smith 
44982bf6240SBarry Smith     Input Parameter:
45082bf6240SBarry Smith .   ts - the timestep context
45182bf6240SBarry Smith 
45282bf6240SBarry Smith     Options Database Key:
45382bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45482bf6240SBarry Smith 
45515091d37SBarry Smith     Level: advanced
45615091d37SBarry Smith 
45782bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
45882bf6240SBarry Smith 
45982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
46082bf6240SBarry Smith @*/
46163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt(TS ts)
46282bf6240SBarry Smith {
463dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS);
46482bf6240SBarry Smith 
46582bf6240SBarry Smith   PetscFunctionBegin;
466*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
46782bf6240SBarry Smith 
468c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
46982bf6240SBarry Smith   if (f) {
47082bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
47182bf6240SBarry Smith   }
47282bf6240SBarry Smith   PetscFunctionReturn(0);
47382bf6240SBarry Smith }
47482bf6240SBarry Smith 
47582bf6240SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
478ac226902SBarry Smith /*@C
47982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
48082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
48182bf6240SBarry Smith 
48215091d37SBarry Smith    Collective on TS
48315091d37SBarry Smith 
48482bf6240SBarry Smith    Input Parameters:
48515091d37SBarry Smith +  ts - timestep context
48682bf6240SBarry Smith .  dt - function to compute timestep
48715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48882bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48982bf6240SBarry Smith 
49015091d37SBarry Smith    Level: intermediate
491fee21e36SBarry Smith 
49282bf6240SBarry Smith    Calling sequence of func:
49387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
49482bf6240SBarry Smith 
49582bf6240SBarry Smith .  newdt - the newly computed timestep
49682bf6240SBarry Smith .  ctx - [optional] timestep context
49782bf6240SBarry Smith 
49882bf6240SBarry Smith    Notes:
49982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
50082bf6240SBarry Smith    during the timestepping process.
50182bf6240SBarry Smith 
50282bf6240SBarry Smith .keywords: timestep, pseudo, set
50382bf6240SBarry Smith 
50482bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
50582bf6240SBarry Smith @*/
50663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
50782bf6240SBarry Smith {
5086849ba73SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *);
50982bf6240SBarry Smith 
51082bf6240SBarry Smith   PetscFunctionBegin;
511*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
51282bf6240SBarry Smith 
513c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
51482bf6240SBarry Smith   if (f) {
51582bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
51682bf6240SBarry Smith   }
51782bf6240SBarry Smith   PetscFunctionReturn(0);
51882bf6240SBarry Smith }
51982bf6240SBarry Smith 
52082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
52182bf6240SBarry Smith 
522a7cc72afSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
523fb2e594dSBarry Smith EXTERN_C_BEGIN
5244a2ae208SSatish Balay #undef __FUNCT__
5254a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
52663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
52782bf6240SBarry Smith {
52882bf6240SBarry Smith   TS_Pseudo *pseudo;
52982bf6240SBarry Smith 
53082bf6240SBarry Smith   PetscFunctionBegin;
53182bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
53282bf6240SBarry Smith   pseudo->verify      = dt;
53382bf6240SBarry Smith   pseudo->verifyctx   = ctx;
53482bf6240SBarry Smith   PetscFunctionReturn(0);
53582bf6240SBarry Smith }
536fb2e594dSBarry Smith EXTERN_C_END
53782bf6240SBarry Smith 
538fb2e594dSBarry Smith EXTERN_C_BEGIN
5394a2ae208SSatish Balay #undef __FUNCT__
5404a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
54163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
54282bf6240SBarry Smith {
5434bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54482bf6240SBarry Smith 
54582bf6240SBarry Smith   PetscFunctionBegin;
54682bf6240SBarry Smith   pseudo->dt_increment = inc;
54782bf6240SBarry Smith   PetscFunctionReturn(0);
54882bf6240SBarry Smith }
549fb2e594dSBarry Smith EXTERN_C_END
55082bf6240SBarry Smith 
551fb2e594dSBarry Smith EXTERN_C_BEGIN
5524a2ae208SSatish Balay #undef __FUNCT__
5534a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
55463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
55582bf6240SBarry Smith {
5564bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55782bf6240SBarry Smith 
55882bf6240SBarry Smith   PetscFunctionBegin;
5594bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
56082bf6240SBarry Smith   PetscFunctionReturn(0);
56182bf6240SBarry Smith }
562fb2e594dSBarry Smith EXTERN_C_END
56382bf6240SBarry Smith 
5646849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
565fb2e594dSBarry Smith EXTERN_C_BEGIN
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
56863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
56982bf6240SBarry Smith {
5704bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
57182bf6240SBarry Smith 
57282bf6240SBarry Smith   PetscFunctionBegin;
57382bf6240SBarry Smith   pseudo->dt      = dt;
57482bf6240SBarry Smith   pseudo->dtctx   = ctx;
57582bf6240SBarry Smith   PetscFunctionReturn(0);
57682bf6240SBarry Smith }
577fb2e594dSBarry Smith EXTERN_C_END
57882bf6240SBarry Smith 
57982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
58010e6a065SJed Brown /*MC
58110e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
58282bf6240SBarry Smith 
58310e6a065SJed Brown   This method solves equations of the form
58410e6a065SJed Brown 
58510e6a065SJed Brown $    F(X,Xdot) = 0
58610e6a065SJed Brown 
58710e6a065SJed Brown   for steady state using the iteration
58810e6a065SJed Brown 
58910e6a065SJed Brown $    [G'] S = -F(X,0)
59010e6a065SJed Brown $    X += S
59110e6a065SJed Brown 
59210e6a065SJed Brown   where
59310e6a065SJed Brown 
59410e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
59510e6a065SJed Brown 
5966f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5976f2d6a7bSJed Brown   state".  See note below.
59810e6a065SJed Brown 
59910e6a065SJed Brown   Options database keys:
60010e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
60110e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
60210e6a065SJed Brown 
60310e6a065SJed Brown   Level: beginner
60410e6a065SJed Brown 
60510e6a065SJed Brown   References:
60610e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
60710e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
60810e6a065SJed Brown 
60910e6a065SJed Brown   Notes:
6106f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6116f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6126f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6136f2d6a7bSJed Brown 
6146f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6156f2d6a7bSJed Brown 
6166f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6176f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6186f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
61910e6a065SJed Brown 
62010e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
62110e6a065SJed Brown 
62210e6a065SJed Brown M*/
623fb2e594dSBarry Smith EXTERN_C_BEGIN
6244a2ae208SSatish Balay #undef __FUNCT__
6254a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
62663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Pseudo(TS ts)
6272d3f70b5SBarry Smith {
6287bf11e45SBarry Smith   TS_Pseudo      *pseudo;
629dfbe8321SBarry Smith   PetscErrorCode ierr;
6302d3f70b5SBarry Smith 
6313a40ed3dSBarry Smith   PetscFunctionBegin;
632000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
633000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6342d3f70b5SBarry Smith 
6352d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
63629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
6372d3f70b5SBarry Smith   }
6388beb423aSHong Zhang   if (!ts->Arhs) {
63929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
6402d3f70b5SBarry Smith   }
641273d9f13SBarry Smith 
642000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
643000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
644000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6457bf11e45SBarry Smith 
6467bf11e45SBarry Smith   /* create the required nonlinear solver context */
6477adad957SLisandro Dalcin   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
64817efb626SBarry Smith   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
6492d3f70b5SBarry Smith 
65038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6517bf11e45SBarry Smith   ts->data = (void*)pseudo;
6522d3f70b5SBarry Smith 
65328aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6544bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
65528aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
65682bf6240SBarry Smith 
657f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
658e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6590c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
660f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
661e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6620c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
663f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
664e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6650c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6660c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6670c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6680c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6693a40ed3dSBarry Smith   PetscFunctionReturn(0);
6702d3f70b5SBarry Smith }
671fb2e594dSBarry Smith EXTERN_C_END
6722d3f70b5SBarry Smith 
6734a2ae208SSatish Balay #undef __FUNCT__
6744a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
67582bf6240SBarry Smith /*@C
67682bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
67782bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
67828aa8177SBarry Smith 
67915091d37SBarry Smith    Collective on TS
68015091d37SBarry Smith 
68128aa8177SBarry Smith    Input Parameters:
68228aa8177SBarry Smith .  ts - the timestep context
68382bf6240SBarry Smith .  dtctx - unused timestep context
68428aa8177SBarry Smith 
68582bf6240SBarry Smith    Output Parameter:
68682bf6240SBarry Smith .  newdt - the timestep to use for the next step
68728aa8177SBarry Smith 
68815091d37SBarry Smith    Level: advanced
68915091d37SBarry Smith 
69082bf6240SBarry Smith .keywords: timestep, pseudo, default
691564e8f4eSLois Curfman McInnes 
69282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
69328aa8177SBarry Smith @*/
69463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
69528aa8177SBarry Smith {
69682bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
69787828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
698dfbe8321SBarry Smith   PetscErrorCode ierr;
69928aa8177SBarry Smith 
7003a40ed3dSBarry Smith   PetscFunctionBegin;
70182bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
70282bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
70382bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
70482bf6240SBarry Smith     /* first time through so compute initial function norm */
70582bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
70682bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
70782bf6240SBarry Smith   }
70882bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
70982bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
710004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
71182bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
71282bf6240SBarry Smith   } else {
71382bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
71482bf6240SBarry Smith   }
71582bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7163a40ed3dSBarry Smith   PetscFunctionReturn(0);
71728aa8177SBarry Smith }
718