xref: /petsc/src/ts/impls/pseudo/posindep.c (revision f1b97656f5355c1550c55f6f45aadd5df55e69a0)
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 
18cdbf8f93SLisandro Dalcin   PetscReal  fnorm_initial,fnorm;                  /* original and current norm of F(u) */
1987828ca2SBarry Smith   PetscReal  fnorm_previous;
2028aa8177SBarry Smith 
21cdbf8f93SLisandro Dalcin   PetscReal  dt_initial;                    /* initial time-step */
2287828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
23ace3abfcSBarry Smith   PetscBool  increment_dt_from_initial_dt;
247bf11e45SBarry Smith } TS_Pseudo;
252d3f70b5SBarry Smith 
262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
272d3f70b5SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
307bf11e45SBarry Smith /*@
317bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
332d3f70b5SBarry Smith 
3415091d37SBarry Smith     Collective on TS
3515091d37SBarry Smith 
367bf11e45SBarry Smith     Input Parameter:
377bf11e45SBarry Smith .   ts - timestep context
387bf11e45SBarry Smith 
397bf11e45SBarry Smith     Output Parameter:
40fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
41fb4a63b6SLois Curfman McInnes 
4215091d37SBarry Smith     Level: advanced
43564e8f4eSLois Curfman McInnes 
44564e8f4eSLois Curfman McInnes     Notes:
45564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
46564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
47564e8f4eSLois Curfman McInnes 
48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
49564e8f4eSLois Curfman McInnes 
50564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
517bf11e45SBarry Smith @*/
527087cfbeSBarry Smith PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
537bf11e45SBarry Smith {
547bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
55dfbe8321SBarry Smith   PetscErrorCode ierr;
567bf11e45SBarry Smith 
573a40ed3dSBarry Smith   PetscFunctionBegin;
58d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
597bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
60d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
613a40ed3dSBarry Smith   PetscFunctionReturn(0);
627bf11e45SBarry Smith }
637bf11e45SBarry Smith 
647bf11e45SBarry Smith 
657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
687bf11e45SBarry Smith /*@C
69639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
707bf11e45SBarry Smith 
7115091d37SBarry Smith    Collective on TS
7215091d37SBarry Smith 
737bf11e45SBarry Smith    Input Parameters:
7415091d37SBarry Smith +  ts - the timestep context
757bf11e45SBarry Smith .  dtctx - unused timestep context
7615091d37SBarry Smith -  update - latest solution vector
777bf11e45SBarry Smith 
78564e8f4eSLois Curfman McInnes    Output Parameters:
7915091d37SBarry Smith +  newdt - the timestep to use for the next step
8015091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
817bf11e45SBarry Smith 
8215091d37SBarry Smith    Level: advanced
83fee21e36SBarry Smith 
84564e8f4eSLois Curfman McInnes    Note:
85564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
86564e8f4eSLois Curfman McInnes    timestep.
87564e8f4eSLois Curfman McInnes 
88564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
89564e8f4eSLois Curfman McInnes 
90564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
917bf11e45SBarry Smith @*/
927087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
937bf11e45SBarry Smith {
943a40ed3dSBarry Smith   PetscFunctionBegin;
95a7cc72afSBarry Smith   *flag = PETSC_TRUE;
963a40ed3dSBarry Smith   PetscFunctionReturn(0);
977bf11e45SBarry Smith }
987bf11e45SBarry Smith 
997bf11e45SBarry Smith 
1004a2ae208SSatish Balay #undef __FUNCT__
1014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1027bf11e45SBarry Smith /*@
103564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1047bf11e45SBarry Smith 
10515091d37SBarry Smith     Collective on TS
10615091d37SBarry Smith 
107fb4a63b6SLois Curfman McInnes     Input Parameters:
10815091d37SBarry Smith +   ts - timestep context
10915091d37SBarry Smith -   update - latest solution vector
1107bf11e45SBarry Smith 
111fb4a63b6SLois Curfman McInnes     Output Parameters:
11215091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11315091d37SBarry Smith -   flag - indicates if current timestep was ok
1147bf11e45SBarry Smith 
11515091d37SBarry Smith     Level: advanced
116fee21e36SBarry Smith 
117564e8f4eSLois Curfman McInnes     Notes:
118564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
119564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
120564e8f4eSLois Curfman McInnes 
121564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
122564e8f4eSLois Curfman McInnes 
123564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1247bf11e45SBarry Smith @*/
1257087cfbeSBarry Smith PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *flag)
1267bf11e45SBarry Smith {
1277bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
128dfbe8321SBarry Smith   PetscErrorCode ierr;
1297bf11e45SBarry Smith 
1303a40ed3dSBarry Smith   PetscFunctionBegin;
131a7cc72afSBarry Smith   if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);}
1327bf11e45SBarry Smith 
1337bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1347bf11e45SBarry Smith 
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1367bf11e45SBarry Smith }
1377bf11e45SBarry Smith 
1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1397bf11e45SBarry Smith 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
142193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1432d3f70b5SBarry Smith {
144277b19d0SLisandro Dalcin   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
145cdbf8f93SLisandro Dalcin   PetscInt       its,lits,reject;
146cdbf8f93SLisandro Dalcin   PetscBool      stepok;
147cdbf8f93SLisandro Dalcin   PetscReal      next_time_step;
148cdbf8f93SLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
149dfbe8321SBarry Smith   PetscErrorCode ierr;
1502d3f70b5SBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
152cdbf8f93SLisandro Dalcin   if (ts->steps == 0) {
153cdbf8f93SLisandro Dalcin     pseudo->dt_initial = ts->time_step;
154cdbf8f93SLisandro Dalcin   }
155193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
156cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
157cdbf8f93SLisandro Dalcin   ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr);
158cdbf8f93SLisandro Dalcin   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
159cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
160193ac0bcSJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr);
161cdbf8f93SLisandro Dalcin     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
162b850b91aSLisandro Dalcin     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1639f954729SBarry Smith     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
164c9780b6fSBarry Smith     ts->nonlinear_its += its; ts->linear_its += lits;
165cdbf8f93SLisandro Dalcin     ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr);
166193ac0bcSJed Brown     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
167cdbf8f93SLisandro Dalcin     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr);
168cdbf8f93SLisandro Dalcin     if (stepok) break;
169cdbf8f93SLisandro Dalcin   }
170*f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
171193ac0bcSJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
172193ac0bcSJed 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);
173cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1747bf11e45SBarry Smith   }
175cdbf8f93SLisandro Dalcin   if (reject >= ts->max_reject) {
176193ac0bcSJed Brown     ts->reason = TS_DIVERGED_STEP_REJECTED;
177cdbf8f93SLisandro Dalcin     ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
178cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
179193ac0bcSJed Brown   }
180cdbf8f93SLisandro Dalcin   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
181193ac0bcSJed Brown   ts->ptime += ts->time_step;
182cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
1832d3f70b5SBarry Smith   ts->steps++;
1843a40ed3dSBarry Smith   PetscFunctionReturn(0);
1852d3f70b5SBarry Smith }
1862d3f70b5SBarry Smith 
1872d3f70b5SBarry Smith /*------------------------------------------------------------*/
1884a2ae208SSatish Balay #undef __FUNCT__
189277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
190277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1912d3f70b5SBarry Smith {
1927bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
193dfbe8321SBarry Smith   PetscErrorCode ierr;
1942d3f70b5SBarry Smith 
1953a40ed3dSBarry Smith   PetscFunctionBegin;
1966bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
1976bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
1986bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
2002d3f70b5SBarry Smith }
2012d3f70b5SBarry Smith 
202277b19d0SLisandro Dalcin #undef __FUNCT__
203277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
204277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
205277b19d0SLisandro Dalcin {
206277b19d0SLisandro Dalcin   PetscErrorCode ierr;
207277b19d0SLisandro Dalcin 
208277b19d0SLisandro Dalcin   PetscFunctionBegin;
209277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
210277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
211335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
212335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr);
213335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr);
214335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
215277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
216277b19d0SLisandro Dalcin }
2172d3f70b5SBarry Smith 
2182d3f70b5SBarry Smith /*------------------------------------------------------------*/
2192d3f70b5SBarry Smith 
2204a2ae208SSatish Balay #undef __FUNCT__
2216f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2226f2d6a7bSJed Brown /*
2236f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2246f2d6a7bSJed Brown */
2256f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2262d3f70b5SBarry Smith {
2276f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
228193ac0bcSJed Brown   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
229193ac0bcSJed Brown   PetscScalar    *xdot;
230dfbe8321SBarry Smith   PetscErrorCode ierr;
231a7cc72afSBarry Smith   PetscInt       i,n;
2322d3f70b5SBarry Smith 
2333a40ed3dSBarry Smith   PetscFunctionBegin;
234193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
235193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2366f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2376f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2382d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2396f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2402d3f70b5SBarry Smith   }
241193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
242193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2436f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2446f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
2462d3f70b5SBarry Smith }
2472d3f70b5SBarry Smith 
2484a2ae208SSatish Balay #undef __FUNCT__
2490f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2506f2d6a7bSJed Brown /*
2516f2d6a7bSJed Brown     The transient residual is
2526f2d6a7bSJed Brown 
2536f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2546f2d6a7bSJed Brown 
2556f2d6a7bSJed Brown     or for ODE,
2566f2d6a7bSJed Brown 
2576f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2586f2d6a7bSJed Brown 
2596f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2606f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2616f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2626f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2636f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2646f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2656f2d6a7bSJed Brown     residual, and then advances to the next time step.
2666f2d6a7bSJed Brown */
2670f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2682d3f70b5SBarry Smith {
2696f2d6a7bSJed Brown   Vec            Xdot;
270dfbe8321SBarry Smith   PetscErrorCode ierr;
2712d3f70b5SBarry Smith 
2723a40ed3dSBarry Smith   PetscFunctionBegin;
2736f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
274193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2756f2d6a7bSJed Brown   PetscFunctionReturn(0);
2766f2d6a7bSJed Brown }
2772d3f70b5SBarry Smith 
2786f2d6a7bSJed Brown #undef __FUNCT__
2790f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2806f2d6a7bSJed Brown /*
2816f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2826f2d6a7bSJed Brown 
2836f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2846f2d6a7bSJed Brown 
2856f2d6a7bSJed Brown     and for ODE:
2866f2d6a7bSJed Brown 
2876f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2886f2d6a7bSJed Brown */
2890f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2906f2d6a7bSJed Brown {
2916f2d6a7bSJed Brown   Vec            Xdot;
2926f2d6a7bSJed Brown   PetscErrorCode ierr;
2936f2d6a7bSJed Brown 
2946f2d6a7bSJed Brown   PetscFunctionBegin;
2956f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
296193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
2982d3f70b5SBarry Smith }
2992d3f70b5SBarry Smith 
3002d3f70b5SBarry Smith 
3014a2ae208SSatish Balay #undef __FUNCT__
3024a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3036849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3042d3f70b5SBarry Smith {
3057bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
306dfbe8321SBarry Smith   PetscErrorCode ierr;
3072d3f70b5SBarry Smith 
3083a40ed3dSBarry Smith   PetscFunctionBegin;
3097bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3107bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3116f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
3132d3f70b5SBarry Smith }
3142d3f70b5SBarry Smith /*------------------------------------------------------------*/
3152d3f70b5SBarry Smith 
3164a2ae208SSatish Balay #undef __FUNCT__
317a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
318649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3192d3f70b5SBarry Smith {
3207bf11e45SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
321dfbe8321SBarry Smith   PetscErrorCode   ierr;
322649052a6SBarry Smith   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3232d3f70b5SBarry Smith 
3243a40ed3dSBarry Smith   PetscFunctionBegin;
325193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
326193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
327193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
328193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
329193ac0bcSJed Brown   }
330649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
331649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
332649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
3342d3f70b5SBarry Smith }
3352d3f70b5SBarry Smith 
3364a2ae208SSatish Balay #undef __FUNCT__
3374a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3386849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3392d3f70b5SBarry Smith {
3404bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
341dfbe8321SBarry Smith   PetscErrorCode  ierr;
342ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
343649052a6SBarry Smith   PetscViewer     viewer;
3442d3f70b5SBarry Smith 
3453a40ed3dSBarry Smith   PetscFunctionBegin;
346b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
347193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3482d3f70b5SBarry Smith     if (flg) {
349649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
350649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
35128aa8177SBarry Smith     }
35290d69ab7SBarry Smith     flg  = PETSC_FALSE;
353acfcf0e5SJed 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);
354ca4b7087SBarry Smith     if (flg) {
355ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
356ca4b7087SBarry Smith     }
357419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
358d52bd9f3SBarry Smith 
359d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
360b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
3622d3f70b5SBarry Smith }
3632d3f70b5SBarry Smith 
3644a2ae208SSatish Balay #undef __FUNCT__
3654a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3666849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3672d3f70b5SBarry Smith {
368d52bd9f3SBarry Smith   PetscErrorCode ierr;
369d52bd9f3SBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
3723a40ed3dSBarry Smith   PetscFunctionReturn(0);
3732d3f70b5SBarry Smith }
3742d3f70b5SBarry Smith 
37582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
378ac226902SBarry Smith /*@C
37982bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
38082bf6240SBarry Smith    last timestep.
38182bf6240SBarry Smith 
3823f9fe445SBarry Smith    Logically Collective on TS
38315091d37SBarry Smith 
38482bf6240SBarry Smith    Input Parameters:
38515091d37SBarry Smith +  ts - timestep context
38682bf6240SBarry Smith .  dt - user-defined function to verify timestep
38715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
38882bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
38982bf6240SBarry Smith 
39015091d37SBarry Smith    Level: advanced
391fee21e36SBarry Smith 
39282bf6240SBarry Smith    Calling sequence of func:
393ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
39482bf6240SBarry Smith 
39582bf6240SBarry Smith .  update - latest solution vector
39682bf6240SBarry Smith .  ctx - [optional] timestep context
39782bf6240SBarry Smith .  newdt - the timestep to use for the next step
39882bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
39982bf6240SBarry Smith 
40082bf6240SBarry Smith    Notes:
40182bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
40282bf6240SBarry Smith    during the timestepping process.
40382bf6240SBarry Smith 
40482bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
40582bf6240SBarry Smith 
40682bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
40782bf6240SBarry Smith @*/
4087087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
40982bf6240SBarry Smith {
4104ac538c5SBarry Smith   PetscErrorCode ierr;
41182bf6240SBarry Smith 
41282bf6240SBarry Smith   PetscFunctionBegin;
4130700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4144ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
41582bf6240SBarry Smith   PetscFunctionReturn(0);
41682bf6240SBarry Smith }
41782bf6240SBarry Smith 
4184a2ae208SSatish Balay #undef __FUNCT__
4194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
42082bf6240SBarry Smith /*@
42182bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
42282bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
42382bf6240SBarry Smith 
4243f9fe445SBarry Smith    Logically Collective on TS
425fee21e36SBarry Smith 
42615091d37SBarry Smith     Input Parameters:
42715091d37SBarry Smith +   ts - the timestep context
42815091d37SBarry Smith -   inc - the scaling factor >= 1.0
42915091d37SBarry Smith 
43082bf6240SBarry Smith     Options Database Key:
43182bf6240SBarry Smith $    -ts_pseudo_increment <increment>
43282bf6240SBarry Smith 
43315091d37SBarry Smith     Level: advanced
43415091d37SBarry Smith 
43582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
43682bf6240SBarry Smith 
43782bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
43882bf6240SBarry Smith @*/
4397087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
44082bf6240SBarry Smith {
4414ac538c5SBarry Smith   PetscErrorCode ierr;
44282bf6240SBarry Smith 
44382bf6240SBarry Smith   PetscFunctionBegin;
4440700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
445c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4464ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
44782bf6240SBarry Smith   PetscFunctionReturn(0);
44882bf6240SBarry Smith }
44982bf6240SBarry Smith 
4504a2ae208SSatish Balay #undef __FUNCT__
4514a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
45282bf6240SBarry Smith /*@
45382bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
45482bf6240SBarry Smith     is computed via the formula
45582bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
45682bf6240SBarry Smith       rather than the default update,
45782bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
45882bf6240SBarry Smith 
4593f9fe445SBarry Smith    Logically Collective on TS
46015091d37SBarry Smith 
46182bf6240SBarry Smith     Input Parameter:
46282bf6240SBarry Smith .   ts - the timestep context
46382bf6240SBarry Smith 
46482bf6240SBarry Smith     Options Database Key:
46582bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
46682bf6240SBarry Smith 
46715091d37SBarry Smith     Level: advanced
46815091d37SBarry Smith 
46982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
47082bf6240SBarry Smith 
47182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
47282bf6240SBarry Smith @*/
4737087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
47482bf6240SBarry Smith {
4754ac538c5SBarry Smith   PetscErrorCode ierr;
47682bf6240SBarry Smith 
47782bf6240SBarry Smith   PetscFunctionBegin;
4780700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4794ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
48082bf6240SBarry Smith   PetscFunctionReturn(0);
48182bf6240SBarry Smith }
48282bf6240SBarry Smith 
48382bf6240SBarry Smith 
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
486ac226902SBarry Smith /*@C
48782bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
48882bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
48982bf6240SBarry Smith 
4903f9fe445SBarry Smith    Logically Collective on TS
49115091d37SBarry Smith 
49282bf6240SBarry Smith    Input Parameters:
49315091d37SBarry Smith +  ts - timestep context
49482bf6240SBarry Smith .  dt - function to compute timestep
49515091d37SBarry Smith -  ctx - [optional] user-defined context for private data
49682bf6240SBarry Smith          required by the function (may be PETSC_NULL)
49782bf6240SBarry Smith 
49815091d37SBarry Smith    Level: intermediate
499fee21e36SBarry Smith 
50082bf6240SBarry Smith    Calling sequence of func:
50187828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
50282bf6240SBarry Smith 
50382bf6240SBarry Smith .  newdt - the newly computed timestep
50482bf6240SBarry Smith .  ctx - [optional] timestep context
50582bf6240SBarry Smith 
50682bf6240SBarry Smith    Notes:
50782bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
50882bf6240SBarry Smith    during the timestepping process.
50982bf6240SBarry Smith 
51082bf6240SBarry Smith .keywords: timestep, pseudo, set
51182bf6240SBarry Smith 
51282bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
51382bf6240SBarry Smith @*/
5147087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
51582bf6240SBarry Smith {
5164ac538c5SBarry Smith   PetscErrorCode ierr;
51782bf6240SBarry Smith 
51882bf6240SBarry Smith   PetscFunctionBegin;
5190700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5204ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
52182bf6240SBarry Smith   PetscFunctionReturn(0);
52282bf6240SBarry Smith }
52382bf6240SBarry Smith 
52482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
52582bf6240SBarry Smith 
526ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
527fb2e594dSBarry Smith EXTERN_C_BEGIN
5284a2ae208SSatish Balay #undef __FUNCT__
5294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5307087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
53182bf6240SBarry Smith {
53282bf6240SBarry Smith   TS_Pseudo *pseudo;
53382bf6240SBarry Smith 
53482bf6240SBarry Smith   PetscFunctionBegin;
53582bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
53682bf6240SBarry Smith   pseudo->verify      = dt;
53782bf6240SBarry Smith   pseudo->verifyctx   = ctx;
53882bf6240SBarry Smith   PetscFunctionReturn(0);
53982bf6240SBarry Smith }
540fb2e594dSBarry Smith EXTERN_C_END
54182bf6240SBarry Smith 
542fb2e594dSBarry Smith EXTERN_C_BEGIN
5434a2ae208SSatish Balay #undef __FUNCT__
5444a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5457087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
54682bf6240SBarry Smith {
5474bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54882bf6240SBarry Smith 
54982bf6240SBarry Smith   PetscFunctionBegin;
55082bf6240SBarry Smith   pseudo->dt_increment = inc;
55182bf6240SBarry Smith   PetscFunctionReturn(0);
55282bf6240SBarry Smith }
553fb2e594dSBarry Smith EXTERN_C_END
55482bf6240SBarry Smith 
555fb2e594dSBarry Smith EXTERN_C_BEGIN
5564a2ae208SSatish Balay #undef __FUNCT__
5574a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
5587087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
55982bf6240SBarry Smith {
5604bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
56182bf6240SBarry Smith 
56282bf6240SBarry Smith   PetscFunctionBegin;
5634bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
56482bf6240SBarry Smith   PetscFunctionReturn(0);
56582bf6240SBarry Smith }
566fb2e594dSBarry Smith EXTERN_C_END
56782bf6240SBarry Smith 
5686849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
569fb2e594dSBarry Smith EXTERN_C_BEGIN
5704a2ae208SSatish Balay #undef __FUNCT__
5714a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
5727087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
57382bf6240SBarry Smith {
5744bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
57582bf6240SBarry Smith 
57682bf6240SBarry Smith   PetscFunctionBegin;
57782bf6240SBarry Smith   pseudo->dt      = dt;
57882bf6240SBarry Smith   pseudo->dtctx   = ctx;
57982bf6240SBarry Smith   PetscFunctionReturn(0);
58082bf6240SBarry Smith }
581fb2e594dSBarry Smith EXTERN_C_END
58282bf6240SBarry Smith 
58382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
58410e6a065SJed Brown /*MC
58510e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
58682bf6240SBarry Smith 
58710e6a065SJed Brown   This method solves equations of the form
58810e6a065SJed Brown 
58910e6a065SJed Brown $    F(X,Xdot) = 0
59010e6a065SJed Brown 
59110e6a065SJed Brown   for steady state using the iteration
59210e6a065SJed Brown 
59310e6a065SJed Brown $    [G'] S = -F(X,0)
59410e6a065SJed Brown $    X += S
59510e6a065SJed Brown 
59610e6a065SJed Brown   where
59710e6a065SJed Brown 
59810e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
59910e6a065SJed Brown 
6006f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6016f2d6a7bSJed Brown   state".  See note below.
60210e6a065SJed Brown 
60310e6a065SJed Brown   Options database keys:
60410e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
60510e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
60610e6a065SJed Brown 
60710e6a065SJed Brown   Level: beginner
60810e6a065SJed Brown 
60910e6a065SJed Brown   References:
61010e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
61110e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
61210e6a065SJed Brown 
61310e6a065SJed Brown   Notes:
6146f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6156f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6166f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6176f2d6a7bSJed Brown 
6186f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6196f2d6a7bSJed Brown 
6206f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6216f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6226f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
62310e6a065SJed Brown 
62410e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
62510e6a065SJed Brown 
62610e6a065SJed Brown M*/
627fb2e594dSBarry Smith EXTERN_C_BEGIN
6284a2ae208SSatish Balay #undef __FUNCT__
6294a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6307087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6312d3f70b5SBarry Smith {
6327bf11e45SBarry Smith   TS_Pseudo      *pseudo;
633dfbe8321SBarry Smith   PetscErrorCode ierr;
634193ac0bcSJed Brown   SNES           snes;
635193ac0bcSJed Brown   const SNESType stype;
6362d3f70b5SBarry Smith 
6373a40ed3dSBarry Smith   PetscFunctionBegin;
638277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
639000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
640000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6412d3f70b5SBarry Smith 
642000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
643000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
644000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6450f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6460f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6477bf11e45SBarry Smith 
648193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
649193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
650193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6512d3f70b5SBarry Smith 
65238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6537bf11e45SBarry Smith   ts->data = (void*)pseudo;
6542d3f70b5SBarry Smith 
65528aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6564bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
65728aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
658193ac0bcSJed Brown   pseudo->fnorm                        = -1;
65982bf6240SBarry Smith 
660f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
661e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6620c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
663f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
664e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6650c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
666f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
667e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6680c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6690c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6700c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6710c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6723a40ed3dSBarry Smith   PetscFunctionReturn(0);
6732d3f70b5SBarry Smith }
674fb2e594dSBarry Smith EXTERN_C_END
6752d3f70b5SBarry Smith 
6764a2ae208SSatish Balay #undef __FUNCT__
6774a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
67882bf6240SBarry Smith /*@C
67982bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
68082bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
68128aa8177SBarry Smith 
68215091d37SBarry Smith    Collective on TS
68315091d37SBarry Smith 
68428aa8177SBarry Smith    Input Parameters:
68528aa8177SBarry Smith .  ts - the timestep context
68682bf6240SBarry Smith .  dtctx - unused timestep context
68728aa8177SBarry Smith 
68882bf6240SBarry Smith    Output Parameter:
68982bf6240SBarry Smith .  newdt - the timestep to use for the next step
69028aa8177SBarry Smith 
69115091d37SBarry Smith    Level: advanced
69215091d37SBarry Smith 
69382bf6240SBarry Smith .keywords: timestep, pseudo, default
694564e8f4eSLois Curfman McInnes 
69582bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
69628aa8177SBarry Smith @*/
6977087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
69828aa8177SBarry Smith {
69982bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
70087828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
701dfbe8321SBarry Smith   PetscErrorCode ierr;
70228aa8177SBarry Smith 
7033a40ed3dSBarry Smith   PetscFunctionBegin;
704bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
705214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
70682bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
707cdbf8f93SLisandro Dalcin   if (pseudo->fnorm_initial == 0.0) {
70882bf6240SBarry Smith     /* first time through so compute initial function norm */
709cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial = pseudo->fnorm;
71082bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
71182bf6240SBarry Smith   }
71282bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
71382bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
714004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
715cdbf8f93SLisandro Dalcin     *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
71682bf6240SBarry Smith   } else {
71782bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
71882bf6240SBarry Smith   }
71982bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7203a40ed3dSBarry Smith   PetscFunctionReturn(0);
72128aa8177SBarry Smith }
722