xref: /petsc/src/ts/impls/pseudo/posindep.c (revision be5899b337ba0cfa5eb720cdae190eefe60949dd)
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;
133cb9d8021SPierre Barbier de Reuille   *flag = PETSC_TRUE;
134*be5899b3SLisandro Dalcin   if(pseudo->verify) {
1357bf11e45SBarry Smith     ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
136cb9d8021SPierre Barbier de Reuille   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
1387bf11e45SBarry Smith }
1397bf11e45SBarry Smith 
1407bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1417bf11e45SBarry Smith 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
144193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts)
1452d3f70b5SBarry Smith {
146277b19d0SLisandro Dalcin   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
147*be5899b3SLisandro Dalcin   PetscInt            nits,lits,reject;
148cdbf8f93SLisandro Dalcin   PetscBool           stepok;
149*be5899b3SLisandro Dalcin   PetscReal           next_time_step = ts->time_step;
150dfbe8321SBarry Smith   PetscErrorCode      ierr;
1512d3f70b5SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
153bbd56ea5SKarl Rupp   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
154193ac0bcSJed Brown   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
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 = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr);
1590298fd71SBarry Smith     ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr);
160*be5899b3SLisandro Dalcin     ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
161b850b91aSLisandro Dalcin     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
162*be5899b3SLisandro Dalcin     ts->snes_its += nits; ts->ksp_its += lits;
1639be3e283SDebojyoti Ghosh     ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr);
164*be5899b3SLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok);CHKERRQ(ierr);
165*be5899b3SLisandro Dalcin     if (!stepok) {next_time_step = ts->time_step; continue;}
166193ac0bcSJed Brown     pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
167cdbf8f93SLisandro Dalcin     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr);
168cdbf8f93SLisandro Dalcin     if (stepok) break;
169cdbf8f93SLisandro Dalcin   }
170*be5899b3SLisandro Dalcin   if (reject >= ts->max_reject) {
171*be5899b3SLisandro Dalcin     ts->reason = TS_DIVERGED_STEP_REJECTED;
172*be5899b3SLisandro Dalcin     ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
173cdbf8f93SLisandro Dalcin     PetscFunctionReturn(0);
1747bf11e45SBarry Smith   }
175*be5899b3SLisandro Dalcin 
176*be5899b3SLisandro Dalcin   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
177*be5899b3SLisandro Dalcin   ts->ptime += ts->time_step;
178*be5899b3SLisandro Dalcin   ts->time_step = next_time_step;
179*be5899b3SLisandro Dalcin 
1803118ae5eSBarry Smith   if (pseudo->fnorm < 0) {
1813118ae5eSBarry Smith     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
1823118ae5eSBarry Smith     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
1833118ae5eSBarry Smith     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
1843118ae5eSBarry Smith   }
1853118ae5eSBarry Smith   if (pseudo->fnorm < pseudo->fatol) {
1863118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
187*be5899b3SLisandro Dalcin     ierr = PetscInfo3(ts,"Step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr);
1883118ae5eSBarry Smith     PetscFunctionReturn(0);
1893118ae5eSBarry Smith   }
1903118ae5eSBarry Smith   if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) {
1913118ae5eSBarry Smith     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
192*be5899b3SLisandro Dalcin     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);
1933118ae5eSBarry Smith     PetscFunctionReturn(0);
1943118ae5eSBarry Smith   }
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1962d3f70b5SBarry Smith }
1972d3f70b5SBarry Smith 
1982d3f70b5SBarry Smith /*------------------------------------------------------------*/
1994a2ae208SSatish Balay #undef __FUNCT__
200277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo"
201277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts)
2022d3f70b5SBarry Smith {
2037bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
204dfbe8321SBarry Smith   PetscErrorCode ierr;
2052d3f70b5SBarry Smith 
2063a40ed3dSBarry Smith   PetscFunctionBegin;
2076bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
2086bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
2096bf464f9SBarry Smith   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
2103a40ed3dSBarry Smith   PetscFunctionReturn(0);
2112d3f70b5SBarry Smith }
2122d3f70b5SBarry Smith 
213277b19d0SLisandro Dalcin #undef __FUNCT__
214277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo"
215277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts)
216277b19d0SLisandro Dalcin {
217277b19d0SLisandro Dalcin   PetscErrorCode ierr;
218277b19d0SLisandro Dalcin 
219277b19d0SLisandro Dalcin   PetscFunctionBegin;
220277b19d0SLisandro Dalcin   ierr = TSReset_Pseudo(ts);CHKERRQ(ierr);
221277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
222bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr);
223bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr);
224bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
225bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr);
226bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr);
227277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
228277b19d0SLisandro Dalcin }
2292d3f70b5SBarry Smith 
2302d3f70b5SBarry Smith /*------------------------------------------------------------*/
2312d3f70b5SBarry Smith 
2324a2ae208SSatish Balay #undef __FUNCT__
2336f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2346f2d6a7bSJed Brown /*
2356f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2366f2d6a7bSJed Brown */
2376f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2382d3f70b5SBarry Smith {
2396f2d6a7bSJed Brown   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
240193ac0bcSJed Brown   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
241193ac0bcSJed Brown   PetscScalar       *xdot;
242dfbe8321SBarry Smith   PetscErrorCode    ierr;
243a7cc72afSBarry Smith   PetscInt          i,n;
2442d3f70b5SBarry Smith 
2453a40ed3dSBarry Smith   PetscFunctionBegin;
246193ac0bcSJed Brown   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
247193ac0bcSJed Brown   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
2486f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2496f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
250bbd56ea5SKarl Rupp   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
251193ac0bcSJed Brown   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
252193ac0bcSJed Brown   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
2536f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2546f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2553a40ed3dSBarry Smith   PetscFunctionReturn(0);
2562d3f70b5SBarry Smith }
2572d3f70b5SBarry Smith 
2584a2ae208SSatish Balay #undef __FUNCT__
2590f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2606f2d6a7bSJed Brown /*
2616f2d6a7bSJed Brown     The transient residual is
2626f2d6a7bSJed Brown 
2636f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2646f2d6a7bSJed Brown 
2656f2d6a7bSJed Brown     or for ODE,
2666f2d6a7bSJed Brown 
2676f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2686f2d6a7bSJed Brown 
2696f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2706f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2716f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2726f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2736f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2746f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2756f2d6a7bSJed Brown     residual, and then advances to the next time step.
2766f2d6a7bSJed Brown */
2770f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2782d3f70b5SBarry Smith {
2796f2d6a7bSJed Brown   Vec            Xdot;
280dfbe8321SBarry Smith   PetscErrorCode ierr;
2812d3f70b5SBarry Smith 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
2836f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
284193ac0bcSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
2856f2d6a7bSJed Brown   PetscFunctionReturn(0);
2866f2d6a7bSJed Brown }
2872d3f70b5SBarry Smith 
2886f2d6a7bSJed Brown #undef __FUNCT__
2890f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2906f2d6a7bSJed Brown /*
2916f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2926f2d6a7bSJed Brown 
2936f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2946f2d6a7bSJed Brown 
2956f2d6a7bSJed Brown     and for ODE:
2966f2d6a7bSJed Brown 
2976f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2986f2d6a7bSJed Brown */
299d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
3006f2d6a7bSJed Brown {
3016f2d6a7bSJed Brown   Vec            Xdot;
3026f2d6a7bSJed Brown   PetscErrorCode ierr;
3036f2d6a7bSJed Brown 
3046f2d6a7bSJed Brown   PetscFunctionBegin;
3056f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
306d1e9a80fSBarry Smith   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr);
3073a40ed3dSBarry Smith   PetscFunctionReturn(0);
3082d3f70b5SBarry Smith }
3092d3f70b5SBarry Smith 
3102d3f70b5SBarry Smith 
3114a2ae208SSatish Balay #undef __FUNCT__
3124a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
3136849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
3142d3f70b5SBarry Smith {
3157bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
316dfbe8321SBarry Smith   PetscErrorCode ierr;
3172d3f70b5SBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3197bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
3207bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
3216f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
3232d3f70b5SBarry Smith }
3242d3f70b5SBarry Smith /*------------------------------------------------------------*/
3252d3f70b5SBarry Smith 
3264a2ae208SSatish Balay #undef __FUNCT__
327a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
328560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3292d3f70b5SBarry Smith {
3307bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
331dfbe8321SBarry Smith   PetscErrorCode ierr;
332ce94432eSBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
3332d3f70b5SBarry Smith 
3343a40ed3dSBarry Smith   PetscFunctionBegin;
335193ac0bcSJed Brown   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
336193ac0bcSJed Brown     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
337193ac0bcSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
338193ac0bcSJed Brown     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
339193ac0bcSJed Brown   }
340649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3417c8652ddSBarry 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);
342649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
3442d3f70b5SBarry Smith }
3452d3f70b5SBarry Smith 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3484416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts)
3492d3f70b5SBarry Smith {
3504bbc92c1SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
351dfbe8321SBarry Smith   PetscErrorCode ierr;
352ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
353649052a6SBarry Smith   PetscViewer    viewer;
3542d3f70b5SBarry Smith 
3553a40ed3dSBarry Smith   PetscFunctionBegin;
356e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr);
357560360afSLisandro Dalcin   ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL);CHKERRQ(ierr);
3582d3f70b5SBarry Smith   if (flg) {
359ce94432eSBarry Smith     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr);
360649052a6SBarry Smith     ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
36128aa8177SBarry Smith   }
362*be5899b3SLisandro Dalcin   flg  = pseudo->increment_dt_from_initial_dt;
3630298fd71SBarry Smith   ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr);
364*be5899b3SLisandro Dalcin   pseudo->increment_dt_from_initial_dt = flg;
36594ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr);
36694ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr);
3673118ae5eSBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr);
3683118ae5eSBarry Smith   ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr);
369b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3703a40ed3dSBarry Smith   PetscFunctionReturn(0);
3712d3f70b5SBarry Smith }
3722d3f70b5SBarry Smith 
3734a2ae208SSatish Balay #undef __FUNCT__
3744a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3756849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3762d3f70b5SBarry Smith {
377d52bd9f3SBarry Smith   PetscErrorCode ierr;
3783118ae5eSBarry Smith   PetscBool      isascii;
379d52bd9f3SBarry Smith 
3803a40ed3dSBarry Smith   PetscFunctionBegin;
3813118ae5eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
3823118ae5eSBarry Smith   if (isascii) {
3833118ae5eSBarry Smith     TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
3843118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Parameters for pseudo timestepping\n");CHKERRQ(ierr);
3853118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr);
3863118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr);
3873118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr);
3883118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr);
3893118ae5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr);
3903118ae5eSBarry Smith   }
391*be5899b3SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
3923a40ed3dSBarry Smith   PetscFunctionReturn(0);
3932d3f70b5SBarry Smith }
3942d3f70b5SBarry Smith 
39582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3964a2ae208SSatish Balay #undef __FUNCT__
3974a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
398ac226902SBarry Smith /*@C
39982bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
40082bf6240SBarry Smith    last timestep.
40182bf6240SBarry Smith 
4023f9fe445SBarry Smith    Logically Collective on TS
40315091d37SBarry Smith 
40482bf6240SBarry Smith    Input Parameters:
40515091d37SBarry Smith +  ts - timestep context
40682bf6240SBarry Smith .  dt - user-defined function to verify timestep
40715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
4080298fd71SBarry Smith          for the timestep verification routine (may be NULL)
40982bf6240SBarry Smith 
41015091d37SBarry Smith    Level: advanced
411fee21e36SBarry Smith 
41282bf6240SBarry Smith    Calling sequence of func:
413ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
41482bf6240SBarry Smith 
41582bf6240SBarry Smith .  update - latest solution vector
41682bf6240SBarry Smith .  ctx - [optional] timestep context
41782bf6240SBarry Smith .  newdt - the timestep to use for the next step
41882bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
41982bf6240SBarry Smith 
42082bf6240SBarry Smith    Notes:
42182bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
42282bf6240SBarry Smith    during the timestepping process.
42382bf6240SBarry Smith 
42482bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
42582bf6240SBarry Smith 
4268d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
42782bf6240SBarry Smith @*/
4287087cfbeSBarry Smith PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
42982bf6240SBarry Smith {
4304ac538c5SBarry Smith   PetscErrorCode ierr;
43182bf6240SBarry Smith 
43282bf6240SBarry Smith   PetscFunctionBegin;
4330700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4344ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr);
43582bf6240SBarry Smith   PetscFunctionReturn(0);
43682bf6240SBarry Smith }
43782bf6240SBarry Smith 
4384a2ae208SSatish Balay #undef __FUNCT__
4394a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
44082bf6240SBarry Smith /*@
44182bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
4428d359177SBarry Smith     dt when using the TSPseudoTimeStepDefault() routine.
44382bf6240SBarry Smith 
4443f9fe445SBarry Smith    Logically Collective on TS
445fee21e36SBarry Smith 
44615091d37SBarry Smith     Input Parameters:
44715091d37SBarry Smith +   ts - the timestep context
44815091d37SBarry Smith -   inc - the scaling factor >= 1.0
44915091d37SBarry Smith 
45082bf6240SBarry Smith     Options Database Key:
451e1bc860dSBarry Smith .    -ts_pseudo_increment <increment>
45282bf6240SBarry Smith 
45315091d37SBarry Smith     Level: advanced
45415091d37SBarry Smith 
45582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
45682bf6240SBarry Smith 
4578d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
45882bf6240SBarry Smith @*/
4597087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
46082bf6240SBarry Smith {
4614ac538c5SBarry Smith   PetscErrorCode ierr;
46282bf6240SBarry Smith 
46382bf6240SBarry Smith   PetscFunctionBegin;
4640700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
465c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4664ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
46782bf6240SBarry Smith   PetscFunctionReturn(0);
46882bf6240SBarry Smith }
46982bf6240SBarry Smith 
4704a2ae208SSatish Balay #undef __FUNCT__
47186534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep"
47286534af1SJed Brown /*@
47386534af1SJed Brown     TSPseudoSetMaxTimeStep - Sets the maximum time step
4748d359177SBarry Smith     when using the TSPseudoTimeStepDefault() routine.
47586534af1SJed Brown 
47686534af1SJed Brown    Logically Collective on TS
47786534af1SJed Brown 
47886534af1SJed Brown     Input Parameters:
47986534af1SJed Brown +   ts - the timestep context
48086534af1SJed Brown -   maxdt - the maximum time step, use a non-positive value to deactivate
48186534af1SJed Brown 
48286534af1SJed Brown     Options Database Key:
483e1bc860dSBarry Smith .    -ts_pseudo_max_dt <increment>
48486534af1SJed Brown 
48586534af1SJed Brown     Level: advanced
48686534af1SJed Brown 
48786534af1SJed Brown .keywords: timestep, pseudo, set
48886534af1SJed Brown 
4898d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
49086534af1SJed Brown @*/
49186534af1SJed Brown PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
49286534af1SJed Brown {
49386534af1SJed Brown   PetscErrorCode ierr;
49486534af1SJed Brown 
49586534af1SJed Brown   PetscFunctionBegin;
49686534af1SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
49786534af1SJed Brown   PetscValidLogicalCollectiveReal(ts,maxdt,2);
49886534af1SJed Brown   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
49986534af1SJed Brown   PetscFunctionReturn(0);
50086534af1SJed Brown }
50186534af1SJed Brown 
50286534af1SJed Brown #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
50482bf6240SBarry Smith /*@
50582bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
50682bf6240SBarry Smith     is computed via the formula
50782bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
50882bf6240SBarry Smith       rather than the default update,
50982bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
51082bf6240SBarry Smith 
5113f9fe445SBarry Smith    Logically Collective on TS
51215091d37SBarry Smith 
51382bf6240SBarry Smith     Input Parameter:
51482bf6240SBarry Smith .   ts - the timestep context
51582bf6240SBarry Smith 
51682bf6240SBarry Smith     Options Database Key:
517e1bc860dSBarry Smith .    -ts_pseudo_increment_dt_from_initial_dt
51882bf6240SBarry Smith 
51915091d37SBarry Smith     Level: advanced
52015091d37SBarry Smith 
52182bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
52282bf6240SBarry Smith 
5238d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
52482bf6240SBarry Smith @*/
5257087cfbeSBarry Smith PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
52682bf6240SBarry Smith {
5274ac538c5SBarry Smith   PetscErrorCode ierr;
52882bf6240SBarry Smith 
52982bf6240SBarry Smith   PetscFunctionBegin;
5300700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5314ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
53282bf6240SBarry Smith   PetscFunctionReturn(0);
53382bf6240SBarry Smith }
53482bf6240SBarry Smith 
53582bf6240SBarry Smith 
5364a2ae208SSatish Balay #undef __FUNCT__
5374a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
538ac226902SBarry Smith /*@C
53982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
54082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
54182bf6240SBarry Smith 
5423f9fe445SBarry Smith    Logically Collective on TS
54315091d37SBarry Smith 
54482bf6240SBarry Smith    Input Parameters:
54515091d37SBarry Smith +  ts - timestep context
54682bf6240SBarry Smith .  dt - function to compute timestep
54715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
5480298fd71SBarry Smith          required by the function (may be NULL)
54982bf6240SBarry Smith 
55015091d37SBarry Smith    Level: intermediate
551fee21e36SBarry Smith 
55282bf6240SBarry Smith    Calling sequence of func:
55387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
55482bf6240SBarry Smith 
55582bf6240SBarry Smith .  newdt - the newly computed timestep
55682bf6240SBarry Smith .  ctx - [optional] timestep context
55782bf6240SBarry Smith 
55882bf6240SBarry Smith    Notes:
55982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
56082bf6240SBarry Smith    during the timestepping process.
5618d359177SBarry Smith    If not set then TSPseudoTimeStepDefault() is automatically used
56282bf6240SBarry Smith 
56382bf6240SBarry Smith .keywords: timestep, pseudo, set
56482bf6240SBarry Smith 
5658d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
56682bf6240SBarry Smith @*/
5677087cfbeSBarry Smith PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
56882bf6240SBarry Smith {
5694ac538c5SBarry Smith   PetscErrorCode ierr;
57082bf6240SBarry Smith 
57182bf6240SBarry Smith   PetscFunctionBegin;
5720700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5734ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr);
57482bf6240SBarry Smith   PetscFunctionReturn(0);
57582bf6240SBarry Smith }
57682bf6240SBarry Smith 
57782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
57882bf6240SBarry Smith 
579ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
582560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
58382bf6240SBarry Smith {
584*be5899b3SLisandro Dalcin   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
58582bf6240SBarry Smith 
58682bf6240SBarry Smith   PetscFunctionBegin;
58782bf6240SBarry Smith   pseudo->verify    = dt;
58882bf6240SBarry Smith   pseudo->verifyctx = ctx;
58982bf6240SBarry Smith   PetscFunctionReturn(0);
59082bf6240SBarry Smith }
59182bf6240SBarry Smith 
5924a2ae208SSatish Balay #undef __FUNCT__
5934a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
594560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
59582bf6240SBarry Smith {
5964bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
59782bf6240SBarry Smith 
59882bf6240SBarry Smith   PetscFunctionBegin;
59982bf6240SBarry Smith   pseudo->dt_increment = inc;
60082bf6240SBarry Smith   PetscFunctionReturn(0);
60182bf6240SBarry Smith }
60282bf6240SBarry Smith 
6034a2ae208SSatish Balay #undef __FUNCT__
60486534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo"
605560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
60686534af1SJed Brown {
60786534af1SJed Brown   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
60886534af1SJed Brown 
60986534af1SJed Brown   PetscFunctionBegin;
61086534af1SJed Brown   pseudo->dt_max = maxdt;
61186534af1SJed Brown   PetscFunctionReturn(0);
61286534af1SJed Brown }
61386534af1SJed Brown 
61486534af1SJed Brown #undef __FUNCT__
6154a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
616560360afSLisandro Dalcin static PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
61782bf6240SBarry Smith {
6184bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61982bf6240SBarry Smith 
62082bf6240SBarry Smith   PetscFunctionBegin;
6214bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
62282bf6240SBarry Smith   PetscFunctionReturn(0);
62382bf6240SBarry Smith }
62482bf6240SBarry Smith 
6256849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
6264a2ae208SSatish Balay #undef __FUNCT__
6274a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
628560360afSLisandro Dalcin static PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
62982bf6240SBarry Smith {
6304bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
63182bf6240SBarry Smith 
63282bf6240SBarry Smith   PetscFunctionBegin;
63382bf6240SBarry Smith   pseudo->dt    = dt;
63482bf6240SBarry Smith   pseudo->dtctx = ctx;
63582bf6240SBarry Smith   PetscFunctionReturn(0);
63682bf6240SBarry Smith }
63782bf6240SBarry Smith 
63882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
63910e6a065SJed Brown /*MC
64010e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
64182bf6240SBarry Smith 
64210e6a065SJed Brown   This method solves equations of the form
64310e6a065SJed Brown 
64410e6a065SJed Brown $    F(X,Xdot) = 0
64510e6a065SJed Brown 
64610e6a065SJed Brown   for steady state using the iteration
64710e6a065SJed Brown 
64810e6a065SJed Brown $    [G'] S = -F(X,0)
64910e6a065SJed Brown $    X += S
65010e6a065SJed Brown 
65110e6a065SJed Brown   where
65210e6a065SJed Brown 
65310e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
65410e6a065SJed Brown 
6556f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
6566f2d6a7bSJed Brown   state".  See note below.
65710e6a065SJed Brown 
65810e6a065SJed Brown   Options database keys:
65910e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
6603118ae5eSBarry Smith .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
6613118ae5eSBarry Smith .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
6623118ae5eSBarry Smith -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
66310e6a065SJed Brown 
66410e6a065SJed Brown   Level: beginner
66510e6a065SJed Brown 
66610e6a065SJed Brown   References:
66796a0c994SBarry Smith +  1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
66896a0c994SBarry Smith -  2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
66910e6a065SJed Brown 
67010e6a065SJed Brown   Notes:
6716f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6726f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6736f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6746f2d6a7bSJed Brown 
6756f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6766f2d6a7bSJed Brown 
6776f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6786f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6796f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
68010e6a065SJed Brown 
68110e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
68210e6a065SJed Brown 
68310e6a065SJed Brown M*/
6844a2ae208SSatish Balay #undef __FUNCT__
6854a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
6868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
6872d3f70b5SBarry Smith {
6887bf11e45SBarry Smith   TS_Pseudo      *pseudo;
689dfbe8321SBarry Smith   PetscErrorCode ierr;
690193ac0bcSJed Brown   SNES           snes;
69119fd82e9SBarry Smith   SNESType       stype;
6922d3f70b5SBarry Smith 
6933a40ed3dSBarry Smith   PetscFunctionBegin;
694277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Pseudo;
695000e7ae3SMatthew Knepley   ts->ops->destroy        = TSDestroy_Pseudo;
696000e7ae3SMatthew Knepley   ts->ops->view           = TSView_Pseudo;
697000e7ae3SMatthew Knepley   ts->ops->setup          = TSSetUp_Pseudo;
698000e7ae3SMatthew Knepley   ts->ops->step           = TSStep_Pseudo;
699000e7ae3SMatthew Knepley   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
7000f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
7010f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
7027bf11e45SBarry Smith 
703193ac0bcSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
704193ac0bcSJed Brown   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
705193ac0bcSJed Brown   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
7062d3f70b5SBarry Smith 
707b00a9115SJed Brown   ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr);
7087bf11e45SBarry Smith   ts->data = (void*)pseudo;
7092d3f70b5SBarry Smith 
710*be5899b3SLisandro Dalcin   pseudo->dt                           = TSPseudoTimeStepDefault;
711*be5899b3SLisandro Dalcin   pseudo->dtctx                        = NULL;
71228aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
7134bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
714193ac0bcSJed Brown   pseudo->fnorm                        = -1;
715*be5899b3SLisandro Dalcin   pseudo->fnorm_initial                = -1;
716*be5899b3SLisandro Dalcin   pseudo->fnorm_previous               = -1;
7173118ae5eSBarry Smith  #if defined(PETSC_USE_REAL_SINGLE)
7183118ae5eSBarry Smith   pseudo->fatol                        = 1.e-25;
7193118ae5eSBarry Smith   pseudo->frtol                        = 1.e-5;
7203118ae5eSBarry Smith #else
7213118ae5eSBarry Smith   pseudo->fatol                        = 1.e-50;
7223118ae5eSBarry Smith   pseudo->frtol                        = 1.e-12;
7233118ae5eSBarry Smith #endif
724bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
725bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
726bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
727bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
728bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
7302d3f70b5SBarry Smith }
7312d3f70b5SBarry Smith 
7324a2ae208SSatish Balay #undef __FUNCT__
7338d359177SBarry Smith #define __FUNCT__ "TSPseudoTimeStepDefault"
73482bf6240SBarry Smith /*@C
7358d359177SBarry Smith    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
73682bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
73728aa8177SBarry Smith 
73815091d37SBarry Smith    Collective on TS
73915091d37SBarry Smith 
74028aa8177SBarry Smith    Input Parameters:
74128aa8177SBarry Smith .  ts - the timestep context
74282bf6240SBarry Smith .  dtctx - unused timestep context
74328aa8177SBarry Smith 
74482bf6240SBarry Smith    Output Parameter:
74582bf6240SBarry Smith .  newdt - the timestep to use for the next step
74628aa8177SBarry Smith 
74715091d37SBarry Smith    Level: advanced
74815091d37SBarry Smith 
74982bf6240SBarry Smith .keywords: timestep, pseudo, default
750564e8f4eSLois Curfman McInnes 
75182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
75228aa8177SBarry Smith @*/
7538d359177SBarry Smith PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
75428aa8177SBarry Smith {
75582bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
756*be5899b3SLisandro Dalcin   PetscReal      inc = pseudo->dt_increment;
757dfbe8321SBarry Smith   PetscErrorCode ierr;
75828aa8177SBarry Smith 
7593a40ed3dSBarry Smith   PetscFunctionBegin;
760bbd7b040SJed Brown   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
761214bc6a2SJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
76282bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
763*be5899b3SLisandro Dalcin   if (pseudo->fnorm_initial < 0) {
76482bf6240SBarry Smith     /* first time through so compute initial function norm */
765cdbf8f93SLisandro Dalcin     pseudo->fnorm_initial  = pseudo->fnorm;
766*be5899b3SLisandro Dalcin     pseudo->fnorm_previous = pseudo->fnorm;
76782bf6240SBarry Smith   }
768bbd56ea5SKarl Rupp   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
769bbd56ea5SKarl Rupp   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
770*be5899b3SLisandro Dalcin   else                                           *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm;
77186534af1SJed Brown   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
77282bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7733a40ed3dSBarry Smith   PetscFunctionReturn(0);
77428aa8177SBarry Smith }
775