xref: /petsc/src/ts/impls/pseudo/posindep.c (revision bbd56ea5790821d2a217d362e8e9710702952333)
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);
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;
152*bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
153193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
154cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
155cdbf8f93SLisandro Dalcin   ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr);
156cdbf8f93SLisandro Dalcin   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
157cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
158b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
159b8123daeSJed Brown     ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr);
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);
1645ef26d82SJed Brown     ts->snes_its += its; ts->ksp_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   }
170f1b97656SJed 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);
21386534af1SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
214335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr);
215335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
216277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
217277b19d0SLisandro Dalcin }
2182d3f70b5SBarry Smith 
2192d3f70b5SBarry Smith /*------------------------------------------------------------*/
2202d3f70b5SBarry Smith 
2214a2ae208SSatish Balay #undef __FUNCT__
2226f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2236f2d6a7bSJed Brown /*
2246f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2256f2d6a7bSJed Brown */
2266f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2272d3f70b5SBarry Smith {
2286f2d6a7bSJed Brown   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
229193ac0bcSJed Brown   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
230193ac0bcSJed Brown   PetscScalar       *xdot;
231dfbe8321SBarry Smith   PetscErrorCode    ierr;
232a7cc72afSBarry Smith   PetscInt          i,n;
2332d3f70b5SBarry Smith 
2343a40ed3dSBarry Smith   PetscFunctionBegin;
235193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
236193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2376f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2386f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
239*bbd56ea5SKarl Rupp   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
240193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
241193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2426f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2436f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2443a40ed3dSBarry Smith   PetscFunctionReturn(0);
2452d3f70b5SBarry Smith }
2462d3f70b5SBarry Smith 
2474a2ae208SSatish Balay #undef __FUNCT__
2480f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2496f2d6a7bSJed Brown /*
2506f2d6a7bSJed Brown     The transient residual is
2516f2d6a7bSJed Brown 
2526f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2536f2d6a7bSJed Brown 
2546f2d6a7bSJed Brown     or for ODE,
2556f2d6a7bSJed Brown 
2566f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2576f2d6a7bSJed Brown 
2586f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2596f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2606f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2616f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2626f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2636f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2646f2d6a7bSJed Brown     residual, and then advances to the next time step.
2656f2d6a7bSJed Brown */
2660f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2672d3f70b5SBarry Smith {
2686f2d6a7bSJed Brown   Vec            Xdot;
269dfbe8321SBarry Smith   PetscErrorCode ierr;
2702d3f70b5SBarry Smith 
2713a40ed3dSBarry Smith   PetscFunctionBegin;
2726f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
273193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2746f2d6a7bSJed Brown   PetscFunctionReturn(0);
2756f2d6a7bSJed Brown }
2762d3f70b5SBarry Smith 
2776f2d6a7bSJed Brown #undef __FUNCT__
2780f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2796f2d6a7bSJed Brown /*
2806f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2816f2d6a7bSJed Brown 
2826f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2836f2d6a7bSJed Brown 
2846f2d6a7bSJed Brown     and for ODE:
2856f2d6a7bSJed Brown 
2866f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2876f2d6a7bSJed Brown */
2880f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2896f2d6a7bSJed Brown {
2906f2d6a7bSJed Brown   Vec            Xdot;
2916f2d6a7bSJed Brown   PetscErrorCode ierr;
2926f2d6a7bSJed Brown 
2936f2d6a7bSJed Brown   PetscFunctionBegin;
2946f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
295193ac0bcSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
2963a40ed3dSBarry Smith   PetscFunctionReturn(0);
2972d3f70b5SBarry Smith }
2982d3f70b5SBarry Smith 
2992d3f70b5SBarry Smith 
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3026849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3032d3f70b5SBarry Smith {
3047bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
305dfbe8321SBarry Smith   PetscErrorCode ierr;
3062d3f70b5SBarry Smith 
3073a40ed3dSBarry Smith   PetscFunctionBegin;
3087bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3097bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3106f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
3122d3f70b5SBarry Smith }
3132d3f70b5SBarry Smith /*------------------------------------------------------------*/
3142d3f70b5SBarry Smith 
3154a2ae208SSatish Balay #undef __FUNCT__
316a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
317649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3182d3f70b5SBarry Smith {
3197bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
320dfbe8321SBarry Smith   PetscErrorCode ierr;
321649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
3222d3f70b5SBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
324193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
325193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
326193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
327193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
328193ac0bcSJed Brown   }
329649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
330649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
331649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3323a40ed3dSBarry Smith   PetscFunctionReturn(0);
3332d3f70b5SBarry Smith }
3342d3f70b5SBarry Smith 
3354a2ae208SSatish Balay #undef __FUNCT__
3364a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3376849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3382d3f70b5SBarry Smith {
3394bbc92c1SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
340dfbe8321SBarry Smith   PetscErrorCode ierr;
341ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
342649052a6SBarry Smith   PetscViewer    viewer;
3432d3f70b5SBarry Smith 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
345b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
346193ac0bcSJed Brown   ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3472d3f70b5SBarry Smith   if (flg) {
348649052a6SBarry Smith     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
349649052a6SBarry Smith     ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
35028aa8177SBarry Smith   }
35190d69ab7SBarry Smith   flg  = PETSC_FALSE;
352acfcf0e5SJed 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);
353ca4b7087SBarry Smith   if (flg) {
354ca4b7087SBarry Smith     ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
355ca4b7087SBarry Smith   }
356419fbf4bSSatish Balay   ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
35786534af1SJed Brown   ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,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__
45186534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep"
45286534af1SJed Brown /*@
45386534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
45486534af1SJed Brown     when using the TSPseudoDefaultTimeStep() routine.
45586534af1SJed Brown 
45686534af1SJed Brown    Logically Collective on TS
45786534af1SJed Brown 
45886534af1SJed Brown     Input Parameters:
45986534af1SJed Brown +   ts - the timestep context
46086534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
46186534af1SJed Brown 
46286534af1SJed Brown     Options Database Key:
46386534af1SJed Brown $    -ts_pseudo_max_dt <increment>
46486534af1SJed Brown 
46586534af1SJed Brown     Level: advanced
46686534af1SJed Brown 
46786534af1SJed Brown .keywords: timestep, pseudo, set
46886534af1SJed Brown 
46986534af1SJed Brown .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
47086534af1SJed Brown @*/
47186534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
47286534af1SJed Brown {
47386534af1SJed Brown   PetscErrorCode ierr;
47486534af1SJed Brown 
47586534af1SJed Brown   PetscFunctionBegin;
47686534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
47786534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
47886534af1SJed Brown   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
47986534af1SJed Brown   PetscFunctionReturn(0);
48086534af1SJed Brown }
48186534af1SJed Brown 
48286534af1SJed Brown #undef __FUNCT__
4834a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
48482bf6240SBarry Smith /*@
48582bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
48682bf6240SBarry Smith     is computed via the formula
48782bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
48882bf6240SBarry Smith       rather than the default update,
48982bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
49082bf6240SBarry Smith 
4913f9fe445SBarry Smith    Logically Collective on TS
49215091d37SBarry Smith 
49382bf6240SBarry Smith     Input Parameter:
49482bf6240SBarry Smith .   ts - the timestep context
49582bf6240SBarry Smith 
49682bf6240SBarry Smith     Options Database Key:
49782bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
49882bf6240SBarry Smith 
49915091d37SBarry Smith     Level: advanced
50015091d37SBarry Smith 
50182bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
50282bf6240SBarry Smith 
50382bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
50482bf6240SBarry Smith @*/
5057087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
50682bf6240SBarry Smith {
5074ac538c5SBarry Smith   PetscErrorCode ierr;
50882bf6240SBarry Smith 
50982bf6240SBarry Smith   PetscFunctionBegin;
5100700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5114ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
51282bf6240SBarry Smith   PetscFunctionReturn(0);
51382bf6240SBarry Smith }
51482bf6240SBarry Smith 
51582bf6240SBarry Smith 
5164a2ae208SSatish Balay #undef __FUNCT__
5174a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
518ac226902SBarry Smith /*@C
51982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
52082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
52182bf6240SBarry Smith 
5223f9fe445SBarry Smith    Logically Collective on TS
52315091d37SBarry Smith 
52482bf6240SBarry Smith    Input Parameters:
52515091d37SBarry Smith +  ts - timestep context
52682bf6240SBarry Smith .  dt - function to compute timestep
52715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
52882bf6240SBarry Smith          required by the function (may be PETSC_NULL)
52982bf6240SBarry Smith 
53015091d37SBarry Smith    Level: intermediate
531fee21e36SBarry Smith 
53282bf6240SBarry Smith    Calling sequence of func:
53387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
53482bf6240SBarry Smith 
53582bf6240SBarry Smith .  newdt - the newly computed timestep
53682bf6240SBarry Smith .  ctx - [optional] timestep context
53782bf6240SBarry Smith 
53882bf6240SBarry Smith    Notes:
53982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
54082bf6240SBarry Smith    during the timestepping process.
54182bf6240SBarry Smith 
54282bf6240SBarry Smith .keywords: timestep, pseudo, set
54382bf6240SBarry Smith 
54482bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
54582bf6240SBarry Smith @*/
5467087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
54782bf6240SBarry Smith {
5484ac538c5SBarry Smith   PetscErrorCode ierr;
54982bf6240SBarry Smith 
55082bf6240SBarry Smith   PetscFunctionBegin;
5510700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5524ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr);
55382bf6240SBarry Smith   PetscFunctionReturn(0);
55482bf6240SBarry Smith }
55582bf6240SBarry Smith 
55682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
55782bf6240SBarry Smith 
558ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
559fb2e594dSBarry Smith EXTERN_C_BEGIN
5604a2ae208SSatish Balay #undef __FUNCT__
5614a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5627087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
56382bf6240SBarry Smith {
56482bf6240SBarry Smith   TS_Pseudo *pseudo;
56582bf6240SBarry Smith 
56682bf6240SBarry Smith   PetscFunctionBegin;
56782bf6240SBarry Smith   pseudo            = (TS_Pseudo*)ts->data;
56882bf6240SBarry Smith   pseudo->verify    = dt;
56982bf6240SBarry Smith   pseudo->verifyctx = ctx;
57082bf6240SBarry Smith   PetscFunctionReturn(0);
57182bf6240SBarry Smith }
572fb2e594dSBarry Smith EXTERN_C_END
57382bf6240SBarry Smith 
574fb2e594dSBarry Smith EXTERN_C_BEGIN
5754a2ae208SSatish Balay #undef __FUNCT__
5764a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
5777087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
57882bf6240SBarry Smith {
5794bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
58082bf6240SBarry Smith 
58182bf6240SBarry Smith   PetscFunctionBegin;
58282bf6240SBarry Smith   pseudo->dt_increment = inc;
58382bf6240SBarry Smith   PetscFunctionReturn(0);
58482bf6240SBarry Smith }
585fb2e594dSBarry Smith EXTERN_C_END
58682bf6240SBarry Smith 
587fb2e594dSBarry Smith EXTERN_C_BEGIN
5884a2ae208SSatish Balay #undef __FUNCT__
58986534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo"
59086534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
59186534af1SJed Brown {
59286534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
59386534af1SJed Brown 
59486534af1SJed Brown   PetscFunctionBegin;
59586534af1SJed Brown   pseudo->dt_max = maxdt;
59686534af1SJed Brown   PetscFunctionReturn(0);
59786534af1SJed Brown }
59886534af1SJed Brown EXTERN_C_END
59986534af1SJed Brown 
60086534af1SJed Brown EXTERN_C_BEGIN
60186534af1SJed Brown #undef __FUNCT__
6024a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
6037087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
60482bf6240SBarry Smith {
6054bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
60682bf6240SBarry Smith 
60782bf6240SBarry Smith   PetscFunctionBegin;
6084bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
60982bf6240SBarry Smith   PetscFunctionReturn(0);
61082bf6240SBarry Smith }
611fb2e594dSBarry Smith EXTERN_C_END
61282bf6240SBarry Smith 
6136849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
614fb2e594dSBarry Smith EXTERN_C_BEGIN
6154a2ae208SSatish Balay #undef __FUNCT__
6164a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
6177087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
61882bf6240SBarry Smith {
6194bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
62082bf6240SBarry Smith 
62182bf6240SBarry Smith   PetscFunctionBegin;
62282bf6240SBarry Smith   pseudo->dt    = dt;
62382bf6240SBarry Smith   pseudo->dtctx = ctx;
62482bf6240SBarry Smith   PetscFunctionReturn(0);
62582bf6240SBarry Smith }
626fb2e594dSBarry Smith EXTERN_C_END
62782bf6240SBarry Smith 
62882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
62910e6a065SJed Brown /*MC
63010e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
63182bf6240SBarry Smith 
63210e6a065SJed Brown   This method solves equations of the form
63310e6a065SJed Brown 
63410e6a065SJed Brown $    F(X,Xdot) = 0
63510e6a065SJed Brown 
63610e6a065SJed Brown   for steady state using the iteration
63710e6a065SJed Brown 
63810e6a065SJed Brown $    [G'] S = -F(X,0)
63910e6a065SJed Brown $    X += S
64010e6a065SJed Brown 
64110e6a065SJed Brown   where
64210e6a065SJed Brown 
64310e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
64410e6a065SJed Brown 
6456f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6466f2d6a7bSJed Brown   state".  See note below.
64710e6a065SJed Brown 
64810e6a065SJed Brown   Options database keys:
64910e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
65010e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
65110e6a065SJed Brown 
65210e6a065SJed Brown   Level: beginner
65310e6a065SJed Brown 
65410e6a065SJed Brown   References:
65510e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
65610e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
65710e6a065SJed Brown 
65810e6a065SJed Brown   Notes:
6596f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6606f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6616f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6626f2d6a7bSJed Brown 
6636f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6646f2d6a7bSJed Brown 
6656f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6666f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6676f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
66810e6a065SJed Brown 
66910e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
67010e6a065SJed Brown 
67110e6a065SJed Brown M*/
672fb2e594dSBarry Smith EXTERN_C_BEGIN
6734a2ae208SSatish Balay #undef __FUNCT__
6744a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6757087cfbeSBarry Smith PetscErrorCode  TSCreate_Pseudo(TS ts)
6762d3f70b5SBarry Smith {
6777bf11e45SBarry Smith   TS_Pseudo      *pseudo;
678dfbe8321SBarry Smith   PetscErrorCode ierr;
679193ac0bcSJed Brown   SNES           snes;
68019fd82e9SBarry Smith   SNESType       stype;
6812d3f70b5SBarry Smith 
6823a40ed3dSBarry Smith   PetscFunctionBegin;
683277b19d0SLisandro Dalcin   ts->ops->reset   = TSReset_Pseudo;
684000e7ae3SMatthew Knepley   ts->ops->destroy = TSDestroy_Pseudo;
685000e7ae3SMatthew Knepley   ts->ops->view    = TSView_Pseudo;
6862d3f70b5SBarry Smith 
687000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
688000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
689000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
6900f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
6910f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
6927bf11e45SBarry Smith 
693193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
694193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
695193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6962d3f70b5SBarry Smith 
69738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6987bf11e45SBarry Smith   ts->data = (void*)pseudo;
6992d3f70b5SBarry Smith 
70028aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7014bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
70228aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
703193ac0bcSJed Brown   pseudo->fnorm                        = -1;
70482bf6240SBarry Smith 
705f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
706e1311b90SBarry Smith                                            "TSPseudoSetVerifyTimeStep_Pseudo",
7070c97f09dSSatish Balay                                            TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
708f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
709e1311b90SBarry Smith                                            "TSPseudoSetTimeStepIncrement_Pseudo",
7100c97f09dSSatish Balay                                            TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
71186534af1SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",
71286534af1SJed Brown                                            "TSPseudoSetMaxTimeStep_Pseudo",
71386534af1SJed Brown                                            TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
714f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
715e1311b90SBarry Smith                                            "TSPseudoIncrementDtFromInitialDt_Pseudo",
7160c97f09dSSatish Balay                                            TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
7170c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
7180c97f09dSSatish Balay                                            "TSPseudoSetTimeStep_Pseudo",
7190c97f09dSSatish Balay                                            TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
7203a40ed3dSBarry Smith   PetscFunctionReturn(0);
7212d3f70b5SBarry Smith }
722fb2e594dSBarry Smith EXTERN_C_END
7232d3f70b5SBarry Smith 
7244a2ae208SSatish Balay #undef __FUNCT__
7254a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
72682bf6240SBarry Smith /*@C
72782bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
72882bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
72928aa8177SBarry Smith 
73015091d37SBarry Smith    Collective on TS
73115091d37SBarry Smith 
73228aa8177SBarry Smith    Input Parameters:
73328aa8177SBarry Smith .  ts - the timestep context
73482bf6240SBarry Smith .  dtctx - unused timestep context
73528aa8177SBarry Smith 
73682bf6240SBarry Smith    Output Parameter:
73782bf6240SBarry Smith .  newdt - the timestep to use for the next step
73828aa8177SBarry Smith 
73915091d37SBarry Smith    Level: advanced
74015091d37SBarry Smith 
74182bf6240SBarry Smith .keywords: timestep, pseudo, default
742564e8f4eSLois Curfman McInnes 
74382bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
74428aa8177SBarry Smith @*/
7457087cfbeSBarry Smith PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal *newdt,void *dtctx)
74628aa8177SBarry Smith {
74782bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
74887828ca2SBarry Smith   PetscReal      inc     = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
749dfbe8321SBarry Smith   PetscErrorCode ierr;
75028aa8177SBarry Smith 
7513a40ed3dSBarry Smith   PetscFunctionBegin;
752bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
753214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
75482bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
755cdbf8f93SLisandro Dalcin   if (pseudo->fnorm_initial == 0.0) {
75682bf6240SBarry Smith     /* first time through so compute initial function norm */
757cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial = pseudo->fnorm;
75882bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
75982bf6240SBarry Smith   }
760*bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
761*bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
762*bbd56ea5SKarl Rupp   else                                           *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
76386534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
76482bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7653a40ed3dSBarry Smith   PetscFunctionReturn(0);
76628aa8177SBarry Smith }
767