xref: /petsc/src/ts/impls/pseudo/posindep.c (revision cb9d8021f64c572d262e10034bbd2fb15469a2b9)
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"
318d359177SBarry Smith /*@C
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 
438d359177SBarry Smith     Level: developer
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 
518d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), 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__
688d359177SBarry Smith #define __FUNCT__ "TSPseudoVerifyTimeStepDefault"
697bf11e45SBarry Smith /*@C
708d359177SBarry Smith    TSPseudoVerifyTimeStepDefault - 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 @*/
938d359177SBarry Smith PetscErrorCode  TSPseudoVerifyTimeStepDefault(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 
1248d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
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;
1327bf11e45SBarry Smith 
133*cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
134*cb9d8021SPierre Barbier de Reuille   ierr = TSFunctionDomainError(ts,ts->ptime,0,&update,flag);CHKERRQ(ierr);
135*cb9d8021SPierre Barbier de Reuille   if(*flag && pseudo->verify) {
1367bf11e45SBarry Smith     ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
137*cb9d8021SPierre Barbier de Reuille   }
1383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1397bf11e45SBarry Smith }
1407bf11e45SBarry Smith 
1417bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1427bf11e45SBarry Smith 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
145193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1462d3f70b5SBarry Smith {
147277b19d0SLisandro Dalcin   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
148cdbf8f93SLisandro Dalcin   PetscInt            its,lits,reject;
149cdbf8f93SLisandro Dalcin   PetscBool           stepok;
150cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
151cdbf8f93SLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
152dfbe8321SBarry Smith   PetscErrorCode      ierr;
1532d3f70b5SBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
155bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
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);
1630298fd71SBarry Smith     ierr = SNESSolve(ts->snes,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);
1679be3e283SDebojyoti Ghosh     ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr);
1685ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
169cdbf8f93SLisandro Dalcin     ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr);
170193ac0bcSJed Brown     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
171cdbf8f93SLisandro Dalcin     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr);
172cdbf8f93SLisandro Dalcin     if (stepok) break;
173cdbf8f93SLisandro Dalcin   }
174f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
175193ac0bcSJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
176193ac0bcSJed 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);
177cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1787bf11e45SBarry Smith   }
179cdbf8f93SLisandro Dalcin   if (reject >= ts->max_reject) {
180193ac0bcSJed Brown     ts->reason = TS_DIVERGED_STEP_REJECTED;
181cdbf8f93SLisandro Dalcin     ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
182cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
183193ac0bcSJed Brown   }
184cdbf8f93SLisandro Dalcin   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
185193ac0bcSJed Brown   ts->ptime += ts->time_step;
186cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
1872d3f70b5SBarry Smith   ts->steps++;
1883a40ed3dSBarry Smith   PetscFunctionReturn(0);
1892d3f70b5SBarry Smith }
1902d3f70b5SBarry Smith 
1912d3f70b5SBarry Smith /*------------------------------------------------------------*/
1924a2ae208SSatish Balay #undef __FUNCT__
193277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
194277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
1952d3f70b5SBarry Smith {
1967bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
197dfbe8321SBarry Smith   PetscErrorCode ierr;
1982d3f70b5SBarry Smith 
1993a40ed3dSBarry Smith   PetscFunctionBegin;
2006bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
2016bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
2026bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
2042d3f70b5SBarry Smith }
2052d3f70b5SBarry Smith 
206277b19d0SLisandro Dalcin #undef __FUNCT__
207277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
208277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
209277b19d0SLisandro Dalcin {
210277b19d0SLisandro Dalcin   PetscErrorCode ierr;
211277b19d0SLisandro Dalcin 
212277b19d0SLisandro Dalcin   PetscFunctionBegin;
213277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
214277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
215bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr);
216bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr);
217bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
218bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr);
219bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr);
220277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
221277b19d0SLisandro Dalcin }
2222d3f70b5SBarry Smith 
2232d3f70b5SBarry Smith /*------------------------------------------------------------*/
2242d3f70b5SBarry Smith 
2254a2ae208SSatish Balay #undef __FUNCT__
2266f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2276f2d6a7bSJed Brown /*
2286f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2296f2d6a7bSJed Brown */
2306f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2312d3f70b5SBarry Smith {
2326f2d6a7bSJed Brown   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
233193ac0bcSJed Brown   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
234193ac0bcSJed Brown   PetscScalar       *xdot;
235dfbe8321SBarry Smith   PetscErrorCode    ierr;
236a7cc72afSBarry Smith   PetscInt          i,n;
2372d3f70b5SBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
239193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
240193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2416f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2426f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
243bbd56ea5SKarl Rupp   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
244193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
245193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2466f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2476f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
2492d3f70b5SBarry Smith }
2502d3f70b5SBarry Smith 
2514a2ae208SSatish Balay #undef __FUNCT__
2520f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2536f2d6a7bSJed Brown /*
2546f2d6a7bSJed Brown     The transient residual is
2556f2d6a7bSJed Brown 
2566f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2576f2d6a7bSJed Brown 
2586f2d6a7bSJed Brown     or for ODE,
2596f2d6a7bSJed Brown 
2606f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2616f2d6a7bSJed Brown 
2626f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2636f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2646f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2656f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2666f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2676f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2686f2d6a7bSJed Brown     residual, and then advances to the next time step.
2696f2d6a7bSJed Brown */
2700f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2712d3f70b5SBarry Smith {
2726f2d6a7bSJed Brown   Vec            Xdot;
273dfbe8321SBarry Smith   PetscErrorCode ierr;
2742d3f70b5SBarry Smith 
2753a40ed3dSBarry Smith   PetscFunctionBegin;
2766f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
277193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2786f2d6a7bSJed Brown   PetscFunctionReturn(0);
2796f2d6a7bSJed Brown }
2802d3f70b5SBarry Smith 
2816f2d6a7bSJed Brown #undef __FUNCT__
2820f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2836f2d6a7bSJed Brown /*
2846f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2856f2d6a7bSJed Brown 
2866f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2876f2d6a7bSJed Brown 
2886f2d6a7bSJed Brown     and for ODE:
2896f2d6a7bSJed Brown 
2906f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2916f2d6a7bSJed Brown */
292d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
2936f2d6a7bSJed Brown {
2946f2d6a7bSJed Brown   Vec            Xdot;
2956f2d6a7bSJed Brown   PetscErrorCode ierr;
2966f2d6a7bSJed Brown 
2976f2d6a7bSJed Brown   PetscFunctionBegin;
2986f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
299d1e9a80fSBarry Smith   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr);
3003a40ed3dSBarry Smith   PetscFunctionReturn(0);
3012d3f70b5SBarry Smith }
3022d3f70b5SBarry Smith 
3032d3f70b5SBarry Smith 
3044a2ae208SSatish Balay #undef __FUNCT__
3054a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3066849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3072d3f70b5SBarry Smith {
3087bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
309dfbe8321SBarry Smith   PetscErrorCode ierr;
3102d3f70b5SBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionBegin;
3127bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3137bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3146f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
3162d3f70b5SBarry Smith }
3172d3f70b5SBarry Smith /*------------------------------------------------------------*/
3182d3f70b5SBarry Smith 
3194a2ae208SSatish Balay #undef __FUNCT__
320a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
321649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3222d3f70b5SBarry Smith {
3237bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
324dfbe8321SBarry Smith   PetscErrorCode ierr;
325ce94432eSBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
3262d3f70b5SBarry Smith 
3273a40ed3dSBarry Smith   PetscFunctionBegin;
328ce94432eSBarry Smith   if (!viewer) {
329ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr);
330ce94432eSBarry Smith   }
331193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
332193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
333193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
334193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
335193ac0bcSJed Brown   }
336649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3377c8652ddSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);CHKERRQ(ierr);
338649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
3402d3f70b5SBarry Smith }
3412d3f70b5SBarry Smith 
3424a2ae208SSatish Balay #undef __FUNCT__
3434a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3448c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptions *PetscOptionsObject,TS ts)
3452d3f70b5SBarry Smith {
3464bbc92c1SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
347dfbe8321SBarry Smith   PetscErrorCode ierr;
348ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
349649052a6SBarry Smith   PetscViewer    viewer;
3502d3f70b5SBarry Smith 
3513a40ed3dSBarry Smith   PetscFunctionBegin;
352e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr);
3530298fd71SBarry Smith   ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);CHKERRQ(ierr);
3542d3f70b5SBarry Smith   if (flg) {
355ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr);
356649052a6SBarry Smith     ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
35728aa8177SBarry Smith   }
35890d69ab7SBarry Smith   flg  = PETSC_FALSE;
3590298fd71SBarry Smith   ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr);
360ca4b7087SBarry Smith   if (flg) {
361ca4b7087SBarry Smith     ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
362ca4b7087SBarry Smith   }
36394ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr);
36494ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr);
365d52bd9f3SBarry Smith 
366d52bd9f3SBarry Smith   ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
367b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3683a40ed3dSBarry Smith   PetscFunctionReturn(0);
3692d3f70b5SBarry Smith }
3702d3f70b5SBarry Smith 
3714a2ae208SSatish Balay #undef __FUNCT__
3724a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3736849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3742d3f70b5SBarry Smith {
375d52bd9f3SBarry Smith   PetscErrorCode ierr;
376d52bd9f3SBarry Smith 
3773a40ed3dSBarry Smith   PetscFunctionBegin;
378d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
3802d3f70b5SBarry Smith }
3812d3f70b5SBarry Smith 
38282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3834a2ae208SSatish Balay #undef __FUNCT__
3844a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
385ac226902SBarry Smith /*@C
38682bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
38782bf6240SBarry Smith    last timestep.
38882bf6240SBarry Smith 
3893f9fe445SBarry Smith    Logically Collective on TS
39015091d37SBarry Smith 
39182bf6240SBarry Smith    Input Parameters:
39215091d37SBarry Smith +  ts - timestep context
39382bf6240SBarry Smith .  dt - user-defined function to verify timestep
39415091d37SBarry Smith -  ctx - [optional] user-defined context for private data
3950298fd71SBarry Smith          for the timestep verification routine (may be NULL)
39682bf6240SBarry Smith 
39715091d37SBarry Smith    Level: advanced
398fee21e36SBarry Smith 
39982bf6240SBarry Smith    Calling sequence of func:
400ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
40182bf6240SBarry Smith 
40282bf6240SBarry Smith .  update - latest solution vector
40382bf6240SBarry Smith .  ctx - [optional] timestep context
40482bf6240SBarry Smith .  newdt - the timestep to use for the next step
40582bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
40682bf6240SBarry Smith 
40782bf6240SBarry Smith    Notes:
40882bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
40982bf6240SBarry Smith    during the timestepping process.
41082bf6240SBarry Smith 
41182bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
41282bf6240SBarry Smith 
4138d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
41482bf6240SBarry Smith @*/
4157087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
41682bf6240SBarry Smith {
4174ac538c5SBarry Smith   PetscErrorCode ierr;
41882bf6240SBarry Smith 
41982bf6240SBarry Smith   PetscFunctionBegin;
4200700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4214ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr);
42282bf6240SBarry Smith   PetscFunctionReturn(0);
42382bf6240SBarry Smith }
42482bf6240SBarry Smith 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
42782bf6240SBarry Smith /*@
42882bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
4298d359177SBarry Smith     dt when using the TSPseudoTimeStepDefault() routine.
43082bf6240SBarry Smith 
4313f9fe445SBarry Smith    Logically Collective on TS
432fee21e36SBarry Smith 
43315091d37SBarry Smith     Input Parameters:
43415091d37SBarry Smith +   ts - the timestep context
43515091d37SBarry Smith -   inc - the scaling factor >= 1.0
43615091d37SBarry Smith 
43782bf6240SBarry Smith     Options Database Key:
438e1bc860dSBarry Smith .    -ts_pseudo_increment <increment>
43982bf6240SBarry Smith 
44015091d37SBarry Smith     Level: advanced
44115091d37SBarry Smith 
44282bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
44382bf6240SBarry Smith 
4448d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
44582bf6240SBarry Smith @*/
4467087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
44782bf6240SBarry Smith {
4484ac538c5SBarry Smith   PetscErrorCode ierr;
44982bf6240SBarry Smith 
45082bf6240SBarry Smith   PetscFunctionBegin;
4510700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
452c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4534ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
45482bf6240SBarry Smith   PetscFunctionReturn(0);
45582bf6240SBarry Smith }
45682bf6240SBarry Smith 
4574a2ae208SSatish Balay #undef __FUNCT__
45886534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep"
45986534af1SJed Brown /*@
46086534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
4618d359177SBarry Smith     when using the TSPseudoTimeStepDefault() routine.
46286534af1SJed Brown 
46386534af1SJed Brown    Logically Collective on TS
46486534af1SJed Brown 
46586534af1SJed Brown     Input Parameters:
46686534af1SJed Brown +   ts - the timestep context
46786534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
46886534af1SJed Brown 
46986534af1SJed Brown     Options Database Key:
470e1bc860dSBarry Smith .    -ts_pseudo_max_dt <increment>
47186534af1SJed Brown 
47286534af1SJed Brown     Level: advanced
47386534af1SJed Brown 
47486534af1SJed Brown .keywords: timestep, pseudo, set
47586534af1SJed Brown 
4768d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
47786534af1SJed Brown @*/
47886534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
47986534af1SJed Brown {
48086534af1SJed Brown   PetscErrorCode ierr;
48186534af1SJed Brown 
48286534af1SJed Brown   PetscFunctionBegin;
48386534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
48486534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
48586534af1SJed Brown   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
48686534af1SJed Brown   PetscFunctionReturn(0);
48786534af1SJed Brown }
48886534af1SJed Brown 
48986534af1SJed Brown #undef __FUNCT__
4904a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
49182bf6240SBarry Smith /*@
49282bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
49382bf6240SBarry Smith     is computed via the formula
49482bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
49582bf6240SBarry Smith       rather than the default update,
49682bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
49782bf6240SBarry Smith 
4983f9fe445SBarry Smith    Logically Collective on TS
49915091d37SBarry Smith 
50082bf6240SBarry Smith     Input Parameter:
50182bf6240SBarry Smith .   ts - the timestep context
50282bf6240SBarry Smith 
50382bf6240SBarry Smith     Options Database Key:
504e1bc860dSBarry Smith .    -ts_pseudo_increment_dt_from_initial_dt
50582bf6240SBarry Smith 
50615091d37SBarry Smith     Level: advanced
50715091d37SBarry Smith 
50882bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
50982bf6240SBarry Smith 
5108d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
51182bf6240SBarry Smith @*/
5127087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
51382bf6240SBarry Smith {
5144ac538c5SBarry Smith   PetscErrorCode ierr;
51582bf6240SBarry Smith 
51682bf6240SBarry Smith   PetscFunctionBegin;
5170700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5184ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
51982bf6240SBarry Smith   PetscFunctionReturn(0);
52082bf6240SBarry Smith }
52182bf6240SBarry Smith 
52282bf6240SBarry Smith 
5234a2ae208SSatish Balay #undef __FUNCT__
5244a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
525ac226902SBarry Smith /*@C
52682bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
52782bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
52882bf6240SBarry Smith 
5293f9fe445SBarry Smith    Logically Collective on TS
53015091d37SBarry Smith 
53182bf6240SBarry Smith    Input Parameters:
53215091d37SBarry Smith +  ts - timestep context
53382bf6240SBarry Smith .  dt - function to compute timestep
53415091d37SBarry Smith -  ctx - [optional] user-defined context for private data
5350298fd71SBarry Smith          required by the function (may be NULL)
53682bf6240SBarry Smith 
53715091d37SBarry Smith    Level: intermediate
538fee21e36SBarry Smith 
53982bf6240SBarry Smith    Calling sequence of func:
54087828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
54182bf6240SBarry Smith 
54282bf6240SBarry Smith .  newdt - the newly computed timestep
54382bf6240SBarry Smith .  ctx - [optional] timestep context
54482bf6240SBarry Smith 
54582bf6240SBarry Smith    Notes:
54682bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
54782bf6240SBarry Smith    during the timestepping process.
5488d359177SBarry Smith    If not set then TSPseudoTimeStepDefault() is automatically used
54982bf6240SBarry Smith 
55082bf6240SBarry Smith .keywords: timestep, pseudo, set
55182bf6240SBarry Smith 
5528d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
55382bf6240SBarry Smith @*/
5547087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
55582bf6240SBarry Smith {
5564ac538c5SBarry Smith   PetscErrorCode ierr;
55782bf6240SBarry Smith 
55882bf6240SBarry Smith   PetscFunctionBegin;
5590700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5604ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr);
56182bf6240SBarry Smith   PetscFunctionReturn(0);
56282bf6240SBarry Smith }
56382bf6240SBarry Smith 
56482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
56582bf6240SBarry Smith 
566ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
5674a2ae208SSatish Balay #undef __FUNCT__
5684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
5697087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
57082bf6240SBarry Smith {
57182bf6240SBarry Smith   TS_Pseudo *pseudo;
57282bf6240SBarry Smith 
57382bf6240SBarry Smith   PetscFunctionBegin;
57482bf6240SBarry Smith   pseudo            = (TS_Pseudo*)ts->data;
57582bf6240SBarry Smith   pseudo->verify    = dt;
57682bf6240SBarry Smith   pseudo->verifyctx = ctx;
57782bf6240SBarry Smith   PetscFunctionReturn(0);
57882bf6240SBarry Smith }
57982bf6240SBarry Smith 
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 }
59082bf6240SBarry Smith 
5914a2ae208SSatish Balay #undef __FUNCT__
59286534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo"
59386534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
59486534af1SJed Brown {
59586534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
59686534af1SJed Brown 
59786534af1SJed Brown   PetscFunctionBegin;
59886534af1SJed Brown   pseudo->dt_max = maxdt;
59986534af1SJed Brown   PetscFunctionReturn(0);
60086534af1SJed Brown }
60186534af1SJed Brown 
60286534af1SJed Brown #undef __FUNCT__
6034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
6047087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
60582bf6240SBarry Smith {
6064bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
60782bf6240SBarry Smith 
60882bf6240SBarry Smith   PetscFunctionBegin;
6094bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
61082bf6240SBarry Smith   PetscFunctionReturn(0);
61182bf6240SBarry Smith }
61282bf6240SBarry Smith 
6136849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
6144a2ae208SSatish Balay #undef __FUNCT__
6154a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
6167087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
61782bf6240SBarry Smith {
6184bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61982bf6240SBarry Smith 
62082bf6240SBarry Smith   PetscFunctionBegin;
62182bf6240SBarry Smith   pseudo->dt    = dt;
62282bf6240SBarry Smith   pseudo->dtctx = ctx;
62382bf6240SBarry Smith   PetscFunctionReturn(0);
62482bf6240SBarry Smith }
62582bf6240SBarry Smith 
62682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
62710e6a065SJed Brown /*MC
62810e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
62982bf6240SBarry Smith 
63010e6a065SJed Brown   This method solves equations of the form
63110e6a065SJed Brown 
63210e6a065SJed Brown $    F(X,Xdot) = 0
63310e6a065SJed Brown 
63410e6a065SJed Brown   for steady state using the iteration
63510e6a065SJed Brown 
63610e6a065SJed Brown $    [G'] S = -F(X,0)
63710e6a065SJed Brown $    X += S
63810e6a065SJed Brown 
63910e6a065SJed Brown   where
64010e6a065SJed Brown 
64110e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
64210e6a065SJed Brown 
6436f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6446f2d6a7bSJed Brown   state".  See note below.
64510e6a065SJed Brown 
64610e6a065SJed Brown   Options database keys:
64710e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
64810e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
64910e6a065SJed Brown 
65010e6a065SJed Brown   Level: beginner
65110e6a065SJed Brown 
65210e6a065SJed Brown   References:
65310e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
65410e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
65510e6a065SJed Brown 
65610e6a065SJed Brown   Notes:
6576f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6586f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6596f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6606f2d6a7bSJed Brown 
6616f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6626f2d6a7bSJed Brown 
6636f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6646f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6656f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
66610e6a065SJed Brown 
66710e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
66810e6a065SJed Brown 
66910e6a065SJed Brown M*/
6704a2ae208SSatish Balay #undef __FUNCT__
6714a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
6732d3f70b5SBarry Smith {
6747bf11e45SBarry Smith   TS_Pseudo      *pseudo;
675dfbe8321SBarry Smith   PetscErrorCode ierr;
676193ac0bcSJed Brown   SNES           snes;
67719fd82e9SBarry Smith   SNESType       stype;
6782d3f70b5SBarry Smith 
6793a40ed3dSBarry Smith   PetscFunctionBegin;
680277b19d0SLisandro Dalcin   ts->ops->reset   = TSReset_Pseudo;
681000e7ae3SMatthew Knepley   ts->ops->destroy = TSDestroy_Pseudo;
682000e7ae3SMatthew Knepley   ts->ops->view    = TSView_Pseudo;
6832d3f70b5SBarry Smith 
684000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
685000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
686000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
6870f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
6880f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
6897bf11e45SBarry Smith 
690193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
691193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
692193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
6932d3f70b5SBarry Smith 
694b00a9115SJed Brown   ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr);
6957bf11e45SBarry Smith   ts->data = (void*)pseudo;
6962d3f70b5SBarry Smith 
69728aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6984bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
6998d359177SBarry Smith   pseudo->dt                           = TSPseudoTimeStepDefault;
700193ac0bcSJed Brown   pseudo->fnorm                        = -1;
70182bf6240SBarry Smith 
702bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
704bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
705bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
706bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
7073a40ed3dSBarry Smith   PetscFunctionReturn(0);
7082d3f70b5SBarry Smith }
7092d3f70b5SBarry Smith 
7104a2ae208SSatish Balay #undef __FUNCT__
7118d359177SBarry Smith #define __FUNCT__ "TSPseudoTimeStepDefault"
71282bf6240SBarry Smith /*@C
7138d359177SBarry Smith    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
71482bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
71528aa8177SBarry Smith 
71615091d37SBarry Smith    Collective on TS
71715091d37SBarry Smith 
71828aa8177SBarry Smith    Input Parameters:
71928aa8177SBarry Smith .  ts - the timestep context
72082bf6240SBarry Smith .  dtctx - unused timestep context
72128aa8177SBarry Smith 
72282bf6240SBarry Smith    Output Parameter:
72382bf6240SBarry Smith .  newdt - the timestep to use for the next step
72428aa8177SBarry Smith 
72515091d37SBarry Smith    Level: advanced
72615091d37SBarry Smith 
72782bf6240SBarry Smith .keywords: timestep, pseudo, default
728564e8f4eSLois Curfman McInnes 
72982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
73028aa8177SBarry Smith @*/
7318d359177SBarry Smith PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
73228aa8177SBarry Smith {
73382bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
73487828ca2SBarry Smith   PetscReal      inc     = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
735dfbe8321SBarry Smith   PetscErrorCode ierr;
73628aa8177SBarry Smith 
7373a40ed3dSBarry Smith   PetscFunctionBegin;
738bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
739214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
74082bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
741cdbf8f93SLisandro Dalcin   if (pseudo->fnorm_initial == 0.0) {
74282bf6240SBarry Smith     /* first time through so compute initial function norm */
743cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial = pseudo->fnorm;
74482bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
74582bf6240SBarry Smith   }
746bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
747bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
748bbd56ea5SKarl Rupp   else                                           *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
74986534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
75082bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7513a40ed3dSBarry Smith   PetscFunctionReturn(0);
75228aa8177SBarry Smith }
753