xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 19fd82e9cecaa743f7d770680f25ee0b5a8ad815)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4b45d2f2cSJed Brown #include <petsc-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 */
2386534af1SJed Brown   PetscReal  dt_max;                        /* maximum time step */
24ace3abfcSBarry Smith   PetscBool  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 @*/
537087cfbeSBarry Smith PetscErrorCode  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 @*/
937087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *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 @*/
1267087cfbeSBarry Smith PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *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"
143193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1442d3f70b5SBarry Smith {
145277b19d0SLisandro Dalcin   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
146cdbf8f93SLisandro Dalcin   PetscInt       its,lits,reject;
147cdbf8f93SLisandro Dalcin   PetscBool      stepok;
148cdbf8f93SLisandro Dalcin   PetscReal      next_time_step;
149cdbf8f93SLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
150dfbe8321SBarry Smith   PetscErrorCode ierr;
1512d3f70b5SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
153cdbf8f93SLisandro Dalcin   if (ts->steps == 0) {
154cdbf8f93SLisandro Dalcin     pseudo->dt_initial = ts->time_step;
155cdbf8f93SLisandro Dalcin   }
156193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
157cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
158cdbf8f93SLisandro Dalcin   ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr);
159cdbf8f93SLisandro Dalcin   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
160cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
161b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
162b8123daeSJed Brown     ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr);
163193ac0bcSJed Brown     ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr);
164cdbf8f93SLisandro Dalcin     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
165b850b91aSLisandro Dalcin     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1669f954729SBarry Smith     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
1675ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
168cdbf8f93SLisandro Dalcin     ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr);
169193ac0bcSJed Brown     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
170cdbf8f93SLisandro Dalcin     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr);
171cdbf8f93SLisandro Dalcin     if (stepok) break;
172cdbf8f93SLisandro Dalcin   }
173f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
174193ac0bcSJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
175193ac0bcSJed 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);
176cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1777bf11e45SBarry Smith   }
178cdbf8f93SLisandro Dalcin   if (reject >= ts->max_reject) {
179193ac0bcSJed Brown     ts->reason = TS_DIVERGED_STEP_REJECTED;
180cdbf8f93SLisandro Dalcin     ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
181cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
182193ac0bcSJed Brown   }
183cdbf8f93SLisandro Dalcin   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
184193ac0bcSJed Brown   ts->ptime += ts->time_step;
185cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
1862d3f70b5SBarry Smith   ts->steps++;
1873a40ed3dSBarry Smith   PetscFunctionReturn(0);
1882d3f70b5SBarry Smith }
1892d3f70b5SBarry Smith 
1902d3f70b5SBarry Smith /*------------------------------------------------------------*/
1914a2ae208SSatish Balay #undef __FUNCT__
192277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
193277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1942d3f70b5SBarry Smith {
1957bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
196dfbe8321SBarry Smith   PetscErrorCode ierr;
1972d3f70b5SBarry Smith 
1983a40ed3dSBarry Smith   PetscFunctionBegin;
1996bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
2006bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
2016bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
2032d3f70b5SBarry Smith }
2042d3f70b5SBarry Smith 
205277b19d0SLisandro Dalcin #undef __FUNCT__
206277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
207277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
208277b19d0SLisandro Dalcin {
209277b19d0SLisandro Dalcin   PetscErrorCode ierr;
210277b19d0SLisandro Dalcin 
211277b19d0SLisandro Dalcin   PetscFunctionBegin;
212277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
213277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
214335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
215335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr);
21686534af1SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
217335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr);
218335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
219277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
220277b19d0SLisandro Dalcin }
2212d3f70b5SBarry Smith 
2222d3f70b5SBarry Smith /*------------------------------------------------------------*/
2232d3f70b5SBarry Smith 
2244a2ae208SSatish Balay #undef __FUNCT__
2256f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2266f2d6a7bSJed Brown /*
2276f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2286f2d6a7bSJed Brown */
2296f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2302d3f70b5SBarry Smith {
2316f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
232193ac0bcSJed Brown   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
233193ac0bcSJed Brown   PetscScalar    *xdot;
234dfbe8321SBarry Smith   PetscErrorCode ierr;
235a7cc72afSBarry Smith   PetscInt       i,n;
2362d3f70b5SBarry Smith 
2373a40ed3dSBarry Smith   PetscFunctionBegin;
238193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
239193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2406f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2416f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2422d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2436f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2442d3f70b5SBarry Smith   }
245193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
246193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2476f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2486f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
2502d3f70b5SBarry Smith }
2512d3f70b5SBarry Smith 
2524a2ae208SSatish Balay #undef __FUNCT__
2530f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2546f2d6a7bSJed Brown /*
2556f2d6a7bSJed Brown     The transient residual is
2566f2d6a7bSJed Brown 
2576f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2586f2d6a7bSJed Brown 
2596f2d6a7bSJed Brown     or for ODE,
2606f2d6a7bSJed Brown 
2616f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2626f2d6a7bSJed Brown 
2636f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2646f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2656f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2666f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2676f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2686f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2696f2d6a7bSJed Brown     residual, and then advances to the next time step.
2706f2d6a7bSJed Brown */
2710f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2722d3f70b5SBarry Smith {
2736f2d6a7bSJed Brown   Vec            Xdot;
274dfbe8321SBarry Smith   PetscErrorCode ierr;
2752d3f70b5SBarry Smith 
2763a40ed3dSBarry Smith   PetscFunctionBegin;
2776f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
278193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2796f2d6a7bSJed Brown   PetscFunctionReturn(0);
2806f2d6a7bSJed Brown }
2812d3f70b5SBarry Smith 
2826f2d6a7bSJed Brown #undef __FUNCT__
2830f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2846f2d6a7bSJed Brown /*
2856f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2866f2d6a7bSJed Brown 
2876f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2886f2d6a7bSJed Brown 
2896f2d6a7bSJed Brown     and for ODE:
2906f2d6a7bSJed Brown 
2916f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2926f2d6a7bSJed Brown */
2930f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2946f2d6a7bSJed Brown {
2956f2d6a7bSJed Brown   Vec            Xdot;
2966f2d6a7bSJed Brown   PetscErrorCode ierr;
2976f2d6a7bSJed Brown 
2986f2d6a7bSJed Brown   PetscFunctionBegin;
2996f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
300193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
3013a40ed3dSBarry Smith   PetscFunctionReturn(0);
3022d3f70b5SBarry Smith }
3032d3f70b5SBarry Smith 
3042d3f70b5SBarry Smith 
3054a2ae208SSatish Balay #undef __FUNCT__
3064a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3076849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3082d3f70b5SBarry Smith {
3097bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
310dfbe8321SBarry Smith   PetscErrorCode ierr;
3112d3f70b5SBarry Smith 
3123a40ed3dSBarry Smith   PetscFunctionBegin;
3137bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3147bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3156f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3163a40ed3dSBarry Smith   PetscFunctionReturn(0);
3172d3f70b5SBarry Smith }
3182d3f70b5SBarry Smith /*------------------------------------------------------------*/
3192d3f70b5SBarry Smith 
3204a2ae208SSatish Balay #undef __FUNCT__
321a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
322649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3232d3f70b5SBarry Smith {
3247bf11e45SBarry Smith   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
325dfbe8321SBarry Smith   PetscErrorCode   ierr;
326649052a6SBarry Smith   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3272d3f70b5SBarry Smith 
3283a40ed3dSBarry Smith   PetscFunctionBegin;
329193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
330193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
331193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
332193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
333193ac0bcSJed Brown   }
334649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
335649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
336649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3373a40ed3dSBarry Smith   PetscFunctionReturn(0);
3382d3f70b5SBarry Smith }
3392d3f70b5SBarry Smith 
3404a2ae208SSatish Balay #undef __FUNCT__
3414a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3426849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3432d3f70b5SBarry Smith {
3444bbc92c1SBarry Smith   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
345dfbe8321SBarry Smith   PetscErrorCode  ierr;
346ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
347649052a6SBarry Smith   PetscViewer     viewer;
3482d3f70b5SBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
350b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
351193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3522d3f70b5SBarry Smith     if (flg) {
353649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
354649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
35528aa8177SBarry Smith     }
35690d69ab7SBarry Smith     flg  = PETSC_FALSE;
357acfcf0e5SJed 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);
358ca4b7087SBarry Smith     if (flg) {
359ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
360ca4b7087SBarry Smith     }
361419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
36286534af1SJed Brown     ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,0);CHKERRQ(ierr);
363d52bd9f3SBarry Smith 
364d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
365b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3663a40ed3dSBarry Smith   PetscFunctionReturn(0);
3672d3f70b5SBarry Smith }
3682d3f70b5SBarry Smith 
3694a2ae208SSatish Balay #undef __FUNCT__
3704a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3716849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3722d3f70b5SBarry Smith {
373d52bd9f3SBarry Smith   PetscErrorCode ierr;
374d52bd9f3SBarry Smith 
3753a40ed3dSBarry Smith   PetscFunctionBegin;
376d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
3773a40ed3dSBarry Smith   PetscFunctionReturn(0);
3782d3f70b5SBarry Smith }
3792d3f70b5SBarry Smith 
38082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3814a2ae208SSatish Balay #undef __FUNCT__
3824a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
383ac226902SBarry Smith /*@C
38482bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
38582bf6240SBarry Smith    last timestep.
38682bf6240SBarry Smith 
3873f9fe445SBarry Smith    Logically Collective on TS
38815091d37SBarry Smith 
38982bf6240SBarry Smith    Input Parameters:
39015091d37SBarry Smith +  ts - timestep context
39182bf6240SBarry Smith .  dt - user-defined function to verify timestep
39215091d37SBarry Smith -  ctx - [optional] user-defined context for private data
39382bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
39482bf6240SBarry Smith 
39515091d37SBarry Smith    Level: advanced
396fee21e36SBarry Smith 
39782bf6240SBarry Smith    Calling sequence of func:
398ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
39982bf6240SBarry Smith 
40082bf6240SBarry Smith .  update - latest solution vector
40182bf6240SBarry Smith .  ctx - [optional] timestep context
40282bf6240SBarry Smith .  newdt - the timestep to use for the next step
40382bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
40482bf6240SBarry Smith 
40582bf6240SBarry Smith    Notes:
40682bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
40782bf6240SBarry Smith    during the timestepping process.
40882bf6240SBarry Smith 
40982bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
41082bf6240SBarry Smith 
41182bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
41282bf6240SBarry Smith @*/
4137087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
41482bf6240SBarry Smith {
4154ac538c5SBarry Smith   PetscErrorCode ierr;
41682bf6240SBarry Smith 
41782bf6240SBarry Smith   PetscFunctionBegin;
4180700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4194ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
42082bf6240SBarry Smith   PetscFunctionReturn(0);
42182bf6240SBarry Smith }
42282bf6240SBarry Smith 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
42582bf6240SBarry Smith /*@
42682bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
42782bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
42882bf6240SBarry Smith 
4293f9fe445SBarry Smith    Logically Collective on TS
430fee21e36SBarry Smith 
43115091d37SBarry Smith     Input Parameters:
43215091d37SBarry Smith +   ts - the timestep context
43315091d37SBarry Smith -   inc - the scaling factor >= 1.0
43415091d37SBarry Smith 
43582bf6240SBarry Smith     Options Database Key:
43682bf6240SBarry Smith $    -ts_pseudo_increment <increment>
43782bf6240SBarry Smith 
43815091d37SBarry Smith     Level: advanced
43915091d37SBarry Smith 
44082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
44182bf6240SBarry Smith 
44282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
44382bf6240SBarry Smith @*/
4447087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
44582bf6240SBarry Smith {
4464ac538c5SBarry Smith   PetscErrorCode ierr;
44782bf6240SBarry Smith 
44882bf6240SBarry Smith   PetscFunctionBegin;
4490700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
450c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4514ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
45282bf6240SBarry Smith   PetscFunctionReturn(0);
45382bf6240SBarry Smith }
45482bf6240SBarry Smith 
4554a2ae208SSatish Balay #undef __FUNCT__
45686534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep"
45786534af1SJed Brown /*@
45886534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
45986534af1SJed Brown     when using the TSPseudoDefaultTimeStep() routine.
46086534af1SJed Brown 
46186534af1SJed Brown    Logically Collective on TS
46286534af1SJed Brown 
46386534af1SJed Brown     Input Parameters:
46486534af1SJed Brown +   ts - the timestep context
46586534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
46686534af1SJed Brown 
46786534af1SJed Brown     Options Database Key:
46886534af1SJed Brown $    -ts_pseudo_max_dt <increment>
46986534af1SJed Brown 
47086534af1SJed Brown     Level: advanced
47186534af1SJed Brown 
47286534af1SJed Brown .keywords: timestep, pseudo, set
47386534af1SJed Brown 
47486534af1SJed Brown .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
47586534af1SJed Brown @*/
47686534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
47786534af1SJed Brown {
47886534af1SJed Brown   PetscErrorCode ierr;
47986534af1SJed Brown 
48086534af1SJed Brown   PetscFunctionBegin;
48186534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
48286534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
48386534af1SJed Brown   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
48486534af1SJed Brown   PetscFunctionReturn(0);
48586534af1SJed Brown }
48686534af1SJed Brown 
48786534af1SJed Brown #undef __FUNCT__
4884a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
48982bf6240SBarry Smith /*@
49082bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
49182bf6240SBarry Smith     is computed via the formula
49282bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
49382bf6240SBarry Smith       rather than the default update,
49482bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
49582bf6240SBarry Smith 
4963f9fe445SBarry Smith    Logically Collective on TS
49715091d37SBarry Smith 
49882bf6240SBarry Smith     Input Parameter:
49982bf6240SBarry Smith .   ts - the timestep context
50082bf6240SBarry Smith 
50182bf6240SBarry Smith     Options Database Key:
50282bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
50382bf6240SBarry Smith 
50415091d37SBarry Smith     Level: advanced
50515091d37SBarry Smith 
50682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
50782bf6240SBarry Smith 
50882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
50982bf6240SBarry Smith @*/
5107087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
51182bf6240SBarry Smith {
5124ac538c5SBarry Smith   PetscErrorCode ierr;
51382bf6240SBarry Smith 
51482bf6240SBarry Smith   PetscFunctionBegin;
5150700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5164ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
51782bf6240SBarry Smith   PetscFunctionReturn(0);
51882bf6240SBarry Smith }
51982bf6240SBarry Smith 
52082bf6240SBarry Smith 
5214a2ae208SSatish Balay #undef __FUNCT__
5224a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
523ac226902SBarry Smith /*@C
52482bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
52582bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
52682bf6240SBarry Smith 
5273f9fe445SBarry Smith    Logically Collective on TS
52815091d37SBarry Smith 
52982bf6240SBarry Smith    Input Parameters:
53015091d37SBarry Smith +  ts - timestep context
53182bf6240SBarry Smith .  dt - function to compute timestep
53215091d37SBarry Smith -  ctx - [optional] user-defined context for private data
53382bf6240SBarry Smith          required by the function (may be PETSC_NULL)
53482bf6240SBarry Smith 
53515091d37SBarry Smith    Level: intermediate
536fee21e36SBarry Smith 
53782bf6240SBarry Smith    Calling sequence of func:
53887828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
53982bf6240SBarry Smith 
54082bf6240SBarry Smith .  newdt - the newly computed timestep
54182bf6240SBarry Smith .  ctx - [optional] timestep context
54282bf6240SBarry Smith 
54382bf6240SBarry Smith    Notes:
54482bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
54582bf6240SBarry Smith    during the timestepping process.
54682bf6240SBarry Smith 
54782bf6240SBarry Smith .keywords: timestep, pseudo, set
54882bf6240SBarry Smith 
54982bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
55082bf6240SBarry Smith @*/
5517087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
55282bf6240SBarry Smith {
5534ac538c5SBarry Smith   PetscErrorCode ierr;
55482bf6240SBarry Smith 
55582bf6240SBarry Smith   PetscFunctionBegin;
5560700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5574ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
55882bf6240SBarry Smith   PetscFunctionReturn(0);
55982bf6240SBarry Smith }
56082bf6240SBarry Smith 
56182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
56282bf6240SBarry Smith 
563ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
564fb2e594dSBarry Smith EXTERN_C_BEGIN
5654a2ae208SSatish Balay #undef __FUNCT__
5664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5677087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
56882bf6240SBarry Smith {
56982bf6240SBarry Smith   TS_Pseudo *pseudo;
57082bf6240SBarry Smith 
57182bf6240SBarry Smith   PetscFunctionBegin;
57282bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
57382bf6240SBarry Smith   pseudo->verify      = dt;
57482bf6240SBarry Smith   pseudo->verifyctx   = ctx;
57582bf6240SBarry Smith   PetscFunctionReturn(0);
57682bf6240SBarry Smith }
577fb2e594dSBarry Smith EXTERN_C_END
57882bf6240SBarry Smith 
579fb2e594dSBarry Smith EXTERN_C_BEGIN
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5827087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
58382bf6240SBarry Smith {
5844bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
58582bf6240SBarry Smith 
58682bf6240SBarry Smith   PetscFunctionBegin;
58782bf6240SBarry Smith   pseudo->dt_increment = inc;
58882bf6240SBarry Smith   PetscFunctionReturn(0);
58982bf6240SBarry Smith }
590fb2e594dSBarry Smith EXTERN_C_END
59182bf6240SBarry Smith 
592fb2e594dSBarry Smith EXTERN_C_BEGIN
5934a2ae208SSatish Balay #undef __FUNCT__
59486534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo"
59586534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
59686534af1SJed Brown {
59786534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
59886534af1SJed Brown 
59986534af1SJed Brown   PetscFunctionBegin;
60086534af1SJed Brown   pseudo->dt_max = maxdt;
60186534af1SJed Brown   PetscFunctionReturn(0);
60286534af1SJed Brown }
60386534af1SJed Brown EXTERN_C_END
60486534af1SJed Brown 
60586534af1SJed Brown EXTERN_C_BEGIN
60686534af1SJed Brown #undef __FUNCT__
6074a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
6087087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
60982bf6240SBarry Smith {
6104bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61182bf6240SBarry Smith 
61282bf6240SBarry Smith   PetscFunctionBegin;
6134bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
61482bf6240SBarry Smith   PetscFunctionReturn(0);
61582bf6240SBarry Smith }
616fb2e594dSBarry Smith EXTERN_C_END
61782bf6240SBarry Smith 
6186849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
619fb2e594dSBarry Smith EXTERN_C_BEGIN
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
6227087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
62382bf6240SBarry Smith {
6244bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
62582bf6240SBarry Smith 
62682bf6240SBarry Smith   PetscFunctionBegin;
62782bf6240SBarry Smith   pseudo->dt      = dt;
62882bf6240SBarry Smith   pseudo->dtctx   = ctx;
62982bf6240SBarry Smith   PetscFunctionReturn(0);
63082bf6240SBarry Smith }
631fb2e594dSBarry Smith EXTERN_C_END
63282bf6240SBarry Smith 
63382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
63410e6a065SJed Brown /*MC
63510e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
63682bf6240SBarry Smith 
63710e6a065SJed Brown   This method solves equations of the form
63810e6a065SJed Brown 
63910e6a065SJed Brown $    F(X,Xdot) = 0
64010e6a065SJed Brown 
64110e6a065SJed Brown   for steady state using the iteration
64210e6a065SJed Brown 
64310e6a065SJed Brown $    [G'] S = -F(X,0)
64410e6a065SJed Brown $    X += S
64510e6a065SJed Brown 
64610e6a065SJed Brown   where
64710e6a065SJed Brown 
64810e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
64910e6a065SJed Brown 
6506f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6516f2d6a7bSJed Brown   state".  See note below.
65210e6a065SJed Brown 
65310e6a065SJed Brown   Options database keys:
65410e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
65510e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
65610e6a065SJed Brown 
65710e6a065SJed Brown   Level: beginner
65810e6a065SJed Brown 
65910e6a065SJed Brown   References:
66010e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
66110e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
66210e6a065SJed Brown 
66310e6a065SJed Brown   Notes:
6646f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6656f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6666f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6676f2d6a7bSJed Brown 
6686f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6696f2d6a7bSJed Brown 
6706f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6716f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6726f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
67310e6a065SJed Brown 
67410e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
67510e6a065SJed Brown 
67610e6a065SJed Brown M*/
677fb2e594dSBarry Smith EXTERN_C_BEGIN
6784a2ae208SSatish Balay #undef __FUNCT__
6794a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6807087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6812d3f70b5SBarry Smith {
6827bf11e45SBarry Smith   TS_Pseudo      *pseudo;
683dfbe8321SBarry Smith   PetscErrorCode ierr;
684193ac0bcSJed Brown   SNES           snes;
685*19fd82e9SBarry Smith   SNESType stype;
6862d3f70b5SBarry Smith 
6873a40ed3dSBarry Smith   PetscFunctionBegin;
688277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Pseudo;
689000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
690000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6912d3f70b5SBarry Smith 
692000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
693000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
694000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6950f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6960f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6977bf11e45SBarry Smith 
698193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
699193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
700193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
7012d3f70b5SBarry Smith 
70238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
7037bf11e45SBarry Smith   ts->data = (void*)pseudo;
7042d3f70b5SBarry Smith 
70528aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7064bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
70728aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
708193ac0bcSJed Brown   pseudo->fnorm                        = -1;
70982bf6240SBarry Smith 
710f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
711e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
7120c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
713f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
714e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
7150c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
71686534af1SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",
71786534af1SJed Brown                     "TSPseudoSetMaxTimeStep_Pseudo",
71886534af1SJed Brown                      TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
719f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
720e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
7210c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
7220c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
7230c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
7240c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
7253a40ed3dSBarry Smith   PetscFunctionReturn(0);
7262d3f70b5SBarry Smith }
727fb2e594dSBarry Smith EXTERN_C_END
7282d3f70b5SBarry Smith 
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
73182bf6240SBarry Smith /*@C
73282bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
73382bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
73428aa8177SBarry Smith 
73515091d37SBarry Smith    Collective on TS
73615091d37SBarry Smith 
73728aa8177SBarry Smith    Input Parameters:
73828aa8177SBarry Smith .  ts - the timestep context
73982bf6240SBarry Smith .  dtctx - unused timestep context
74028aa8177SBarry Smith 
74182bf6240SBarry Smith    Output Parameter:
74282bf6240SBarry Smith .  newdt - the timestep to use for the next step
74328aa8177SBarry Smith 
74415091d37SBarry Smith    Level: advanced
74515091d37SBarry Smith 
74682bf6240SBarry Smith .keywords: timestep, pseudo, default
747564e8f4eSLois Curfman McInnes 
74882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
74928aa8177SBarry Smith @*/
7507087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
75128aa8177SBarry Smith {
75282bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
75387828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
754dfbe8321SBarry Smith   PetscErrorCode ierr;
75528aa8177SBarry Smith 
7563a40ed3dSBarry Smith   PetscFunctionBegin;
757bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
758214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
75982bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
760cdbf8f93SLisandro Dalcin   if (pseudo->fnorm_initial == 0.0) {
76182bf6240SBarry Smith     /* first time through so compute initial function norm */
762cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial = pseudo->fnorm;
76382bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
76482bf6240SBarry Smith   }
76582bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
76682bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
767004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
768cdbf8f93SLisandro Dalcin     *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
76982bf6240SBarry Smith   } else {
77082bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
77182bf6240SBarry Smith   }
77286534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
77382bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7743a40ed3dSBarry Smith   PetscFunctionReturn(0);
77528aa8177SBarry Smith }
776