xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 560360af1e0197128c3ad6271bbfa2e76ad345d4)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4af0996ceSBarry Smith #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;
253118ae5eSBarry Smith   PetscReal fatol,frtol;
267bf11e45SBarry Smith } TS_Pseudo;
272d3f70b5SBarry Smith 
282d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
292d3f70b5SBarry Smith 
304a2ae208SSatish Balay #undef __FUNCT__
314a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
328d359177SBarry Smith /*@C
337bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
34564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
352d3f70b5SBarry Smith 
3615091d37SBarry Smith     Collective on TS
3715091d37SBarry Smith 
387bf11e45SBarry Smith     Input Parameter:
397bf11e45SBarry Smith .   ts - timestep context
407bf11e45SBarry Smith 
417bf11e45SBarry Smith     Output Parameter:
42fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
43fb4a63b6SLois Curfman McInnes 
448d359177SBarry Smith     Level: developer
45564e8f4eSLois Curfman McInnes 
46564e8f4eSLois Curfman McInnes     Notes:
47564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
48564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
49564e8f4eSLois Curfman McInnes 
50fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
51564e8f4eSLois Curfman McInnes 
528d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep()
537bf11e45SBarry Smith @*/
547087cfbeSBarry Smith PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
557bf11e45SBarry Smith {
567bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
57dfbe8321SBarry Smith   PetscErrorCode ierr;
587bf11e45SBarry Smith 
593a40ed3dSBarry Smith   PetscFunctionBegin;
60d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
617bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
62d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
633a40ed3dSBarry Smith   PetscFunctionReturn(0);
647bf11e45SBarry Smith }
657bf11e45SBarry Smith 
667bf11e45SBarry Smith 
677bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
684a2ae208SSatish Balay #undef __FUNCT__
698d359177SBarry Smith #define __FUNCT__ "TSPseudoVerifyTimeStepDefault"
707bf11e45SBarry Smith /*@C
718d359177SBarry Smith    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
727bf11e45SBarry Smith 
7315091d37SBarry Smith    Collective on TS
7415091d37SBarry Smith 
757bf11e45SBarry Smith    Input Parameters:
7615091d37SBarry Smith +  ts - the timestep context
777bf11e45SBarry Smith .  dtctx - unused timestep context
7815091d37SBarry Smith -  update - latest solution vector
797bf11e45SBarry Smith 
80564e8f4eSLois Curfman McInnes    Output Parameters:
8115091d37SBarry Smith +  newdt - the timestep to use for the next step
8215091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
837bf11e45SBarry Smith 
8415091d37SBarry Smith    Level: advanced
85fee21e36SBarry Smith 
86564e8f4eSLois Curfman McInnes    Note:
87564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
88564e8f4eSLois Curfman McInnes    timestep.
89564e8f4eSLois Curfman McInnes 
90564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
91564e8f4eSLois Curfman McInnes 
92564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
937bf11e45SBarry Smith @*/
948d359177SBarry Smith PetscErrorCode  TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
957bf11e45SBarry Smith {
963a40ed3dSBarry Smith   PetscFunctionBegin;
97a7cc72afSBarry Smith   *flag = PETSC_TRUE;
983a40ed3dSBarry Smith   PetscFunctionReturn(0);
997bf11e45SBarry Smith }
1007bf11e45SBarry Smith 
1017bf11e45SBarry Smith 
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1047bf11e45SBarry Smith /*@
105564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1067bf11e45SBarry Smith 
10715091d37SBarry Smith     Collective on TS
10815091d37SBarry Smith 
109fb4a63b6SLois Curfman McInnes     Input Parameters:
11015091d37SBarry Smith +   ts - timestep context
11115091d37SBarry Smith -   update - latest solution vector
1127bf11e45SBarry Smith 
113fb4a63b6SLois Curfman McInnes     Output Parameters:
11415091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11515091d37SBarry Smith -   flag - indicates if current timestep was ok
1167bf11e45SBarry Smith 
11715091d37SBarry Smith     Level: advanced
118fee21e36SBarry Smith 
119564e8f4eSLois Curfman McInnes     Notes:
120564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
121564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
122564e8f4eSLois Curfman McInnes 
123564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
124564e8f4eSLois Curfman McInnes 
1258d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
1267bf11e45SBarry Smith @*/
1277087cfbeSBarry Smith PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *flag)
1287bf11e45SBarry Smith {
1297bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
130dfbe8321SBarry Smith   PetscErrorCode ierr;
1317bf11e45SBarry Smith 
1323a40ed3dSBarry Smith   PetscFunctionBegin;
1337bf11e45SBarry Smith 
134cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
135d183316bSPierre Barbier de Reuille   ierr = TSFunctionDomainError(ts,ts->ptime,update,flag);CHKERRQ(ierr);
136cb9d8021SPierre Barbier de Reuille   if(*flag && pseudo->verify) {
1377bf11e45SBarry Smith     ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
138cb9d8021SPierre Barbier de Reuille   }
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1407bf11e45SBarry Smith }
1417bf11e45SBarry Smith 
1427bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1437bf11e45SBarry Smith 
1444a2ae208SSatish Balay #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
146193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1472d3f70b5SBarry Smith {
148277b19d0SLisandro Dalcin   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
149cdbf8f93SLisandro Dalcin   PetscInt            its,lits,reject;
150cdbf8f93SLisandro Dalcin   PetscBool           stepok;
151cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
152cdbf8f93SLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
153dfbe8321SBarry Smith   PetscErrorCode      ierr;
1542d3f70b5SBarry Smith 
1553a40ed3dSBarry Smith   PetscFunctionBegin;
156bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
157193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
158cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
159cdbf8f93SLisandro Dalcin   ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr);
160cdbf8f93SLisandro Dalcin   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
161cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
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   }
1793118ae5eSBarry Smith   if (pseudo->fnorm < 0) {
1803118ae5eSBarry Smith     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
1813118ae5eSBarry Smith     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
1823118ae5eSBarry Smith     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
1833118ae5eSBarry Smith   }
1843118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1853118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
1863118ae5eSBarry Smith     ierr = PetscInfo3(ts,"step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr);
1873118ae5eSBarry Smith     PetscFunctionReturn(0);
1883118ae5eSBarry Smith   }
1893118ae5eSBarry Smith   if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) {
1903118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
1913118ae5eSBarry Smith     ierr = PetscInfo4(ts,"step=%D, converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,pseudo->fnorm,pseudo->fnorm_initial,pseudo->fatol);CHKERRQ(ierr);
1923118ae5eSBarry Smith     PetscFunctionReturn(0);
1933118ae5eSBarry Smith   }
194cdbf8f93SLisandro Dalcin   if (reject >= ts->max_reject) {
195193ac0bcSJed Brown     ts->reason = TS_DIVERGED_STEP_REJECTED;
196cdbf8f93SLisandro Dalcin     ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
197cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
198193ac0bcSJed Brown   }
199cdbf8f93SLisandro Dalcin   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
200193ac0bcSJed Brown   ts->ptime += ts->time_step;
201cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
2022d3f70b5SBarry Smith   ts->steps++;
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
2042d3f70b5SBarry Smith }
2052d3f70b5SBarry Smith 
2062d3f70b5SBarry Smith /*------------------------------------------------------------*/
2074a2ae208SSatish Balay #undef __FUNCT__
208277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
209277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
2102d3f70b5SBarry Smith {
2117bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
212dfbe8321SBarry Smith   PetscErrorCode ierr;
2132d3f70b5SBarry Smith 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
2156bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
2166bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
2176bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
2183a40ed3dSBarry Smith   PetscFunctionReturn(0);
2192d3f70b5SBarry Smith }
2202d3f70b5SBarry Smith 
221277b19d0SLisandro Dalcin #undef __FUNCT__
222277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
223277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
224277b19d0SLisandro Dalcin {
225277b19d0SLisandro Dalcin   PetscErrorCode ierr;
226277b19d0SLisandro Dalcin 
227277b19d0SLisandro Dalcin   PetscFunctionBegin;
228277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
229277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
230bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr);
231bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr);
232bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
233bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr);
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr);
235277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
236277b19d0SLisandro Dalcin }
2372d3f70b5SBarry Smith 
2382d3f70b5SBarry Smith /*------------------------------------------------------------*/
2392d3f70b5SBarry Smith 
2404a2ae208SSatish Balay #undef __FUNCT__
2416f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2426f2d6a7bSJed Brown /*
2436f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2446f2d6a7bSJed Brown */
2456f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2462d3f70b5SBarry Smith {
2476f2d6a7bSJed Brown   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
248193ac0bcSJed Brown   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
249193ac0bcSJed Brown   PetscScalar       *xdot;
250dfbe8321SBarry Smith   PetscErrorCode    ierr;
251a7cc72afSBarry Smith   PetscInt          i,n;
2522d3f70b5SBarry Smith 
2533a40ed3dSBarry Smith   PetscFunctionBegin;
254193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
255193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2566f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2576f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
258bbd56ea5SKarl Rupp   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
259193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
260193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2616f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2626f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2633a40ed3dSBarry Smith   PetscFunctionReturn(0);
2642d3f70b5SBarry Smith }
2652d3f70b5SBarry Smith 
2664a2ae208SSatish Balay #undef __FUNCT__
2670f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2686f2d6a7bSJed Brown /*
2696f2d6a7bSJed Brown     The transient residual is
2706f2d6a7bSJed Brown 
2716f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2726f2d6a7bSJed Brown 
2736f2d6a7bSJed Brown     or for ODE,
2746f2d6a7bSJed Brown 
2756f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2766f2d6a7bSJed Brown 
2776f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2786f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2796f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2806f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2816f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2826f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2836f2d6a7bSJed Brown     residual, and then advances to the next time step.
2846f2d6a7bSJed Brown */
2850f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2862d3f70b5SBarry Smith {
2876f2d6a7bSJed Brown   Vec            Xdot;
288dfbe8321SBarry Smith   PetscErrorCode ierr;
2892d3f70b5SBarry Smith 
2903a40ed3dSBarry Smith   PetscFunctionBegin;
2916f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
292193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2936f2d6a7bSJed Brown   PetscFunctionReturn(0);
2946f2d6a7bSJed Brown }
2952d3f70b5SBarry Smith 
2966f2d6a7bSJed Brown #undef __FUNCT__
2970f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2986f2d6a7bSJed Brown /*
2996f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
3006f2d6a7bSJed Brown 
3016f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
3026f2d6a7bSJed Brown 
3036f2d6a7bSJed Brown     and for ODE:
3046f2d6a7bSJed Brown 
3056f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
3066f2d6a7bSJed Brown */
307d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
3086f2d6a7bSJed Brown {
3096f2d6a7bSJed Brown   Vec            Xdot;
3106f2d6a7bSJed Brown   PetscErrorCode ierr;
3116f2d6a7bSJed Brown 
3126f2d6a7bSJed Brown   PetscFunctionBegin;
3136f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
314d1e9a80fSBarry Smith   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
3162d3f70b5SBarry Smith }
3172d3f70b5SBarry Smith 
3182d3f70b5SBarry Smith 
3194a2ae208SSatish Balay #undef __FUNCT__
3204a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3216849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3222d3f70b5SBarry Smith {
3237bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
324dfbe8321SBarry Smith   PetscErrorCode ierr;
3252d3f70b5SBarry Smith 
3263a40ed3dSBarry Smith   PetscFunctionBegin;
3277bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3287bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3296f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3303a40ed3dSBarry Smith   PetscFunctionReturn(0);
3312d3f70b5SBarry Smith }
3322d3f70b5SBarry Smith /*------------------------------------------------------------*/
3332d3f70b5SBarry Smith 
3344a2ae208SSatish Balay #undef __FUNCT__
335a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
336*560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3372d3f70b5SBarry Smith {
3387bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
339dfbe8321SBarry Smith   PetscErrorCode ierr;
340ce94432eSBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
3412d3f70b5SBarry Smith 
3423a40ed3dSBarry Smith   PetscFunctionBegin;
343193ac0bcSJed Brown   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
344193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
345193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
346193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
347193ac0bcSJed Brown   }
348649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3497c8652ddSBarry 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);
350649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3522d3f70b5SBarry Smith }
3532d3f70b5SBarry Smith 
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3564416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts)
3572d3f70b5SBarry Smith {
3584bbc92c1SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
359dfbe8321SBarry Smith   PetscErrorCode ierr;
360ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
361649052a6SBarry Smith   PetscViewer    viewer;
3622d3f70b5SBarry Smith 
3633a40ed3dSBarry Smith   PetscFunctionBegin;
364e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr);
365*560360afSLisandro Dalcin   ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL);CHKERRQ(ierr);
3662d3f70b5SBarry Smith   if (flg) {
367ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr);
368649052a6SBarry Smith     ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
36928aa8177SBarry Smith   }
37090d69ab7SBarry Smith   flg  = PETSC_FALSE;
3710298fd71SBarry Smith   ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr);
372ca4b7087SBarry Smith   if (flg) {
373ca4b7087SBarry Smith     ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
374ca4b7087SBarry Smith   }
37594ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr);
37694ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr);
3773118ae5eSBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr);
3783118ae5eSBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr);
379b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
3812d3f70b5SBarry Smith }
3822d3f70b5SBarry Smith 
3834a2ae208SSatish Balay #undef __FUNCT__
3844a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3856849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3862d3f70b5SBarry Smith {
387d52bd9f3SBarry Smith   PetscErrorCode ierr;
3883118ae5eSBarry Smith   PetscBool      isascii;
389d52bd9f3SBarry Smith 
3903a40ed3dSBarry Smith   PetscFunctionBegin;
3913118ae5eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
3923118ae5eSBarry Smith   if (isascii) {
3933118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3943118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Parameters for pseudo timestepping\n");CHKERRQ(ierr);
3953118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr);
3963118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr);
3973118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr);
3983118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr);
3993118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr);
4003118ae5eSBarry Smith   }
401d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
4032d3f70b5SBarry Smith }
4042d3f70b5SBarry Smith 
40582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
4064a2ae208SSatish Balay #undef __FUNCT__
4074a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
408ac226902SBarry Smith /*@C
40982bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
41082bf6240SBarry Smith    last timestep.
41182bf6240SBarry Smith 
4123f9fe445SBarry Smith    Logically Collective on TS
41315091d37SBarry Smith 
41482bf6240SBarry Smith    Input Parameters:
41515091d37SBarry Smith +  ts - timestep context
41682bf6240SBarry Smith .  dt - user-defined function to verify timestep
41715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
4180298fd71SBarry Smith          for the timestep verification routine (may be NULL)
41982bf6240SBarry Smith 
42015091d37SBarry Smith    Level: advanced
421fee21e36SBarry Smith 
42282bf6240SBarry Smith    Calling sequence of func:
423ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
42482bf6240SBarry Smith 
42582bf6240SBarry Smith .  update - latest solution vector
42682bf6240SBarry Smith .  ctx - [optional] timestep context
42782bf6240SBarry Smith .  newdt - the timestep to use for the next step
42882bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
42982bf6240SBarry Smith 
43082bf6240SBarry Smith    Notes:
43182bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
43282bf6240SBarry Smith    during the timestepping process.
43382bf6240SBarry Smith 
43482bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
43582bf6240SBarry Smith 
4368d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
43782bf6240SBarry Smith @*/
4387087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
43982bf6240SBarry Smith {
4404ac538c5SBarry Smith   PetscErrorCode ierr;
44182bf6240SBarry Smith 
44282bf6240SBarry Smith   PetscFunctionBegin;
4430700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4444ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr);
44582bf6240SBarry Smith   PetscFunctionReturn(0);
44682bf6240SBarry Smith }
44782bf6240SBarry Smith 
4484a2ae208SSatish Balay #undef __FUNCT__
4494a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
45082bf6240SBarry Smith /*@
45182bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
4528d359177SBarry Smith     dt when using the TSPseudoTimeStepDefault() routine.
45382bf6240SBarry Smith 
4543f9fe445SBarry Smith    Logically Collective on TS
455fee21e36SBarry Smith 
45615091d37SBarry Smith     Input Parameters:
45715091d37SBarry Smith +   ts - the timestep context
45815091d37SBarry Smith -   inc - the scaling factor >= 1.0
45915091d37SBarry Smith 
46082bf6240SBarry Smith     Options Database Key:
461e1bc860dSBarry Smith .    -ts_pseudo_increment <increment>
46282bf6240SBarry Smith 
46315091d37SBarry Smith     Level: advanced
46415091d37SBarry Smith 
46582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
46682bf6240SBarry Smith 
4678d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
46882bf6240SBarry Smith @*/
4697087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
47082bf6240SBarry Smith {
4714ac538c5SBarry Smith   PetscErrorCode ierr;
47282bf6240SBarry Smith 
47382bf6240SBarry Smith   PetscFunctionBegin;
4740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
475c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4764ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
47782bf6240SBarry Smith   PetscFunctionReturn(0);
47882bf6240SBarry Smith }
47982bf6240SBarry Smith 
4804a2ae208SSatish Balay #undef __FUNCT__
48186534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep"
48286534af1SJed Brown /*@
48386534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
4848d359177SBarry Smith     when using the TSPseudoTimeStepDefault() routine.
48586534af1SJed Brown 
48686534af1SJed Brown    Logically Collective on TS
48786534af1SJed Brown 
48886534af1SJed Brown     Input Parameters:
48986534af1SJed Brown +   ts - the timestep context
49086534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
49186534af1SJed Brown 
49286534af1SJed Brown     Options Database Key:
493e1bc860dSBarry Smith .    -ts_pseudo_max_dt <increment>
49486534af1SJed Brown 
49586534af1SJed Brown     Level: advanced
49686534af1SJed Brown 
49786534af1SJed Brown .keywords: timestep, pseudo, set
49886534af1SJed Brown 
4998d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
50086534af1SJed Brown @*/
50186534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
50286534af1SJed Brown {
50386534af1SJed Brown   PetscErrorCode ierr;
50486534af1SJed Brown 
50586534af1SJed Brown   PetscFunctionBegin;
50686534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
50786534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
50886534af1SJed Brown   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
50986534af1SJed Brown   PetscFunctionReturn(0);
51086534af1SJed Brown }
51186534af1SJed Brown 
51286534af1SJed Brown #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
51482bf6240SBarry Smith /*@
51582bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
51682bf6240SBarry Smith     is computed via the formula
51782bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
51882bf6240SBarry Smith       rather than the default update,
51982bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
52082bf6240SBarry Smith 
5213f9fe445SBarry Smith    Logically Collective on TS
52215091d37SBarry Smith 
52382bf6240SBarry Smith     Input Parameter:
52482bf6240SBarry Smith .   ts - the timestep context
52582bf6240SBarry Smith 
52682bf6240SBarry Smith     Options Database Key:
527e1bc860dSBarry Smith .    -ts_pseudo_increment_dt_from_initial_dt
52882bf6240SBarry Smith 
52915091d37SBarry Smith     Level: advanced
53015091d37SBarry Smith 
53182bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
53282bf6240SBarry Smith 
5338d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
53482bf6240SBarry Smith @*/
5357087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
53682bf6240SBarry Smith {
5374ac538c5SBarry Smith   PetscErrorCode ierr;
53882bf6240SBarry Smith 
53982bf6240SBarry Smith   PetscFunctionBegin;
5400700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5414ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
54282bf6240SBarry Smith   PetscFunctionReturn(0);
54382bf6240SBarry Smith }
54482bf6240SBarry Smith 
54582bf6240SBarry Smith 
5464a2ae208SSatish Balay #undef __FUNCT__
5474a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
548ac226902SBarry Smith /*@C
54982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
55082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
55182bf6240SBarry Smith 
5523f9fe445SBarry Smith    Logically Collective on TS
55315091d37SBarry Smith 
55482bf6240SBarry Smith    Input Parameters:
55515091d37SBarry Smith +  ts - timestep context
55682bf6240SBarry Smith .  dt - function to compute timestep
55715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
5580298fd71SBarry Smith          required by the function (may be NULL)
55982bf6240SBarry Smith 
56015091d37SBarry Smith    Level: intermediate
561fee21e36SBarry Smith 
56282bf6240SBarry Smith    Calling sequence of func:
56387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
56482bf6240SBarry Smith 
56582bf6240SBarry Smith .  newdt - the newly computed timestep
56682bf6240SBarry Smith .  ctx - [optional] timestep context
56782bf6240SBarry Smith 
56882bf6240SBarry Smith    Notes:
56982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
57082bf6240SBarry Smith    during the timestepping process.
5718d359177SBarry Smith    If not set then TSPseudoTimeStepDefault() is automatically used
57282bf6240SBarry Smith 
57382bf6240SBarry Smith .keywords: timestep, pseudo, set
57482bf6240SBarry Smith 
5758d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
57682bf6240SBarry Smith @*/
5777087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
57882bf6240SBarry Smith {
5794ac538c5SBarry Smith   PetscErrorCode ierr;
58082bf6240SBarry Smith 
58182bf6240SBarry Smith   PetscFunctionBegin;
5820700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5834ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr);
58482bf6240SBarry Smith   PetscFunctionReturn(0);
58582bf6240SBarry Smith }
58682bf6240SBarry Smith 
58782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
58882bf6240SBarry Smith 
589ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
5904a2ae208SSatish Balay #undef __FUNCT__
5914a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
592*560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
59382bf6240SBarry Smith {
59482bf6240SBarry Smith   TS_Pseudo *pseudo;
59582bf6240SBarry Smith 
59682bf6240SBarry Smith   PetscFunctionBegin;
59782bf6240SBarry Smith   pseudo            = (TS_Pseudo*)ts->data;
59882bf6240SBarry Smith   pseudo->verify    = dt;
59982bf6240SBarry Smith   pseudo->verifyctx = ctx;
60082bf6240SBarry Smith   PetscFunctionReturn(0);
60182bf6240SBarry Smith }
60282bf6240SBarry Smith 
6034a2ae208SSatish Balay #undef __FUNCT__
6044a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
605*560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
60682bf6240SBarry Smith {
6074bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
60882bf6240SBarry Smith 
60982bf6240SBarry Smith   PetscFunctionBegin;
61082bf6240SBarry Smith   pseudo->dt_increment = inc;
61182bf6240SBarry Smith   PetscFunctionReturn(0);
61282bf6240SBarry Smith }
61382bf6240SBarry Smith 
6144a2ae208SSatish Balay #undef __FUNCT__
61586534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo"
616*560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
61786534af1SJed Brown {
61886534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61986534af1SJed Brown 
62086534af1SJed Brown   PetscFunctionBegin;
62186534af1SJed Brown   pseudo->dt_max = maxdt;
62286534af1SJed Brown   PetscFunctionReturn(0);
62386534af1SJed Brown }
62486534af1SJed Brown 
62586534af1SJed Brown #undef __FUNCT__
6264a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
627*560360afSLisandro Dalcin static PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
62882bf6240SBarry Smith {
6294bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
63082bf6240SBarry Smith 
63182bf6240SBarry Smith   PetscFunctionBegin;
6324bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
63382bf6240SBarry Smith   PetscFunctionReturn(0);
63482bf6240SBarry Smith }
63582bf6240SBarry Smith 
6366849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
6374a2ae208SSatish Balay #undef __FUNCT__
6384a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
639*560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
64082bf6240SBarry Smith {
6414bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
64282bf6240SBarry Smith 
64382bf6240SBarry Smith   PetscFunctionBegin;
64482bf6240SBarry Smith   pseudo->dt    = dt;
64582bf6240SBarry Smith   pseudo->dtctx = ctx;
64682bf6240SBarry Smith   PetscFunctionReturn(0);
64782bf6240SBarry Smith }
64882bf6240SBarry Smith 
64982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
65010e6a065SJed Brown /*MC
65110e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
65282bf6240SBarry Smith 
65310e6a065SJed Brown   This method solves equations of the form
65410e6a065SJed Brown 
65510e6a065SJed Brown $    F(X,Xdot) = 0
65610e6a065SJed Brown 
65710e6a065SJed Brown   for steady state using the iteration
65810e6a065SJed Brown 
65910e6a065SJed Brown $    [G'] S = -F(X,0)
66010e6a065SJed Brown $    X += S
66110e6a065SJed Brown 
66210e6a065SJed Brown   where
66310e6a065SJed Brown 
66410e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
66510e6a065SJed Brown 
6666f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6676f2d6a7bSJed Brown   state".  See note below.
66810e6a065SJed Brown 
66910e6a065SJed Brown   Options database keys:
67010e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
6713118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
6723118ae5eSBarry Smith .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
6733118ae5eSBarry Smith -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
67410e6a065SJed Brown 
67510e6a065SJed Brown   Level: beginner
67610e6a065SJed Brown 
67710e6a065SJed Brown   References:
67896a0c994SBarry Smith +  1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
67996a0c994SBarry Smith -  2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
68010e6a065SJed Brown 
68110e6a065SJed Brown   Notes:
6826f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6836f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6846f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6856f2d6a7bSJed Brown 
6866f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6876f2d6a7bSJed Brown 
6886f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6896f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6906f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
69110e6a065SJed Brown 
69210e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
69310e6a065SJed Brown 
69410e6a065SJed Brown M*/
6954a2ae208SSatish Balay #undef __FUNCT__
6964a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
6982d3f70b5SBarry Smith {
6997bf11e45SBarry Smith   TS_Pseudo      *pseudo;
700dfbe8321SBarry Smith   PetscErrorCode ierr;
701193ac0bcSJed Brown   SNES           snes;
70219fd82e9SBarry Smith   SNESType       stype;
7032d3f70b5SBarry Smith 
7043a40ed3dSBarry Smith   PetscFunctionBegin;
705277b19d0SLisandro Dalcin   ts->ops->reset   = TSReset_Pseudo;
706000e7ae3SMatthew Knepley   ts->ops->destroy = TSDestroy_Pseudo;
707000e7ae3SMatthew Knepley   ts->ops->view    = TSView_Pseudo;
7082d3f70b5SBarry Smith 
709000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
710000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
711000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
7120f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
7130f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
7147bf11e45SBarry Smith 
715193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
716193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
717193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
7182d3f70b5SBarry Smith 
719b00a9115SJed Brown   ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr);
7207bf11e45SBarry Smith   ts->data = (void*)pseudo;
7212d3f70b5SBarry Smith 
72228aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7234bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
7248d359177SBarry Smith   pseudo->dt                           = TSPseudoTimeStepDefault;
725193ac0bcSJed Brown   pseudo->fnorm                        = -1;
7263118ae5eSBarry Smith  #if defined(PETSC_USE_REAL_SINGLE)
7273118ae5eSBarry Smith   pseudo->fatol                        = 1.e-25;
7283118ae5eSBarry Smith   pseudo->frtol                        = 1.e-5;
7293118ae5eSBarry Smith #else
7303118ae5eSBarry Smith   pseudo->fatol                        = 1.e-50;
7313118ae5eSBarry Smith   pseudo->frtol                        = 1.e-12;
7323118ae5eSBarry Smith #endif
733bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
734bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
735bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
736bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
737bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
7383a40ed3dSBarry Smith   PetscFunctionReturn(0);
7392d3f70b5SBarry Smith }
7402d3f70b5SBarry Smith 
7414a2ae208SSatish Balay #undef __FUNCT__
7428d359177SBarry Smith #define __FUNCT__ "TSPseudoTimeStepDefault"
74382bf6240SBarry Smith /*@C
7448d359177SBarry Smith    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
74582bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
74628aa8177SBarry Smith 
74715091d37SBarry Smith    Collective on TS
74815091d37SBarry Smith 
74928aa8177SBarry Smith    Input Parameters:
75028aa8177SBarry Smith .  ts - the timestep context
75182bf6240SBarry Smith .  dtctx - unused timestep context
75228aa8177SBarry Smith 
75382bf6240SBarry Smith    Output Parameter:
75482bf6240SBarry Smith .  newdt - the timestep to use for the next step
75528aa8177SBarry Smith 
75615091d37SBarry Smith    Level: advanced
75715091d37SBarry Smith 
75882bf6240SBarry Smith .keywords: timestep, pseudo, default
759564e8f4eSLois Curfman McInnes 
76082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
76128aa8177SBarry Smith @*/
7628d359177SBarry Smith PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
76328aa8177SBarry Smith {
76482bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
76587828ca2SBarry Smith   PetscReal      inc     = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
766dfbe8321SBarry Smith   PetscErrorCode ierr;
76728aa8177SBarry Smith 
7683a40ed3dSBarry Smith   PetscFunctionBegin;
769bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
770214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
77182bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
772cdbf8f93SLisandro Dalcin   if (pseudo->fnorm_initial == 0.0) {
77382bf6240SBarry Smith     /* first time through so compute initial function norm */
774cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial = pseudo->fnorm;
77582bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
77682bf6240SBarry Smith   }
777bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
778bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
779bbd56ea5SKarl Rupp   else                                           *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
78086534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
78182bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7823a40ed3dSBarry Smith   PetscFunctionReturn(0);
78328aa8177SBarry Smith }
784