xref: /petsc/src/ts/impls/pseudo/posindep.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
32d3f70b5SBarry Smith /*
4fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
52d3f70b5SBarry Smith */
67c4f633dSBarry Smith #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
72d3f70b5SBarry Smith 
82d3f70b5SBarry Smith typedef struct {
92d3f70b5SBarry Smith   Vec  update;      /* work vector where new solution is formed */
102d3f70b5SBarry Smith   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
116f2d6a7bSJed Brown   Vec  xdot;        /* work vector for time derivative of state */
122d3f70b5SBarry Smith 
132d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
142d3f70b5SBarry Smith 
156849ba73SBarry Smith   PetscErrorCode (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
162d3f70b5SBarry Smith   void           *dtctx;
17ace3abfcSBarry Smith   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscBool *); /* verify previous timestep and related context */
187bf11e45SBarry Smith   void           *verifyctx;
192d3f70b5SBarry Smith 
2087828ca2SBarry Smith   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
2187828ca2SBarry Smith   PetscReal  fnorm_previous;
2228aa8177SBarry Smith 
2387828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
24ace3abfcSBarry Smith   PetscBool  increment_dt_from_initial_dt;
257bf11e45SBarry Smith } TS_Pseudo;
262d3f70b5SBarry Smith 
272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
282d3f70b5SBarry Smith 
294a2ae208SSatish Balay #undef __FUNCT__
304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
317bf11e45SBarry Smith /*@
327bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
33564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
342d3f70b5SBarry Smith 
3515091d37SBarry Smith     Collective on TS
3615091d37SBarry Smith 
377bf11e45SBarry Smith     Input Parameter:
387bf11e45SBarry Smith .   ts - timestep context
397bf11e45SBarry Smith 
407bf11e45SBarry Smith     Output Parameter:
41fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
42fb4a63b6SLois Curfman McInnes 
4315091d37SBarry Smith     Level: advanced
44564e8f4eSLois Curfman McInnes 
45564e8f4eSLois Curfman McInnes     Notes:
46564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
47564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
48564e8f4eSLois Curfman McInnes 
49fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
50564e8f4eSLois Curfman McInnes 
51564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
527bf11e45SBarry Smith @*/
5363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
547bf11e45SBarry Smith {
557bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
56dfbe8321SBarry Smith   PetscErrorCode ierr;
577bf11e45SBarry Smith 
583a40ed3dSBarry Smith   PetscFunctionBegin;
59d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
607bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
61d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
623a40ed3dSBarry Smith   PetscFunctionReturn(0);
637bf11e45SBarry Smith }
647bf11e45SBarry Smith 
657bf11e45SBarry Smith 
667bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
697bf11e45SBarry Smith /*@C
70639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
717bf11e45SBarry Smith 
7215091d37SBarry Smith    Collective on TS
7315091d37SBarry Smith 
747bf11e45SBarry Smith    Input Parameters:
7515091d37SBarry Smith +  ts - the timestep context
767bf11e45SBarry Smith .  dtctx - unused timestep context
7715091d37SBarry Smith -  update - latest solution vector
787bf11e45SBarry Smith 
79564e8f4eSLois Curfman McInnes    Output Parameters:
8015091d37SBarry Smith +  newdt - the timestep to use for the next step
8115091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
827bf11e45SBarry Smith 
8315091d37SBarry Smith    Level: advanced
84fee21e36SBarry Smith 
85564e8f4eSLois Curfman McInnes    Note:
86564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
87564e8f4eSLois Curfman McInnes    timestep.
88564e8f4eSLois Curfman McInnes 
89564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
90564e8f4eSLois Curfman McInnes 
91564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
927bf11e45SBarry Smith @*/
93ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
947bf11e45SBarry Smith {
953a40ed3dSBarry Smith   PetscFunctionBegin;
96a7cc72afSBarry Smith   *flag = PETSC_TRUE;
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
987bf11e45SBarry Smith }
997bf11e45SBarry Smith 
1007bf11e45SBarry Smith 
1014a2ae208SSatish Balay #undef __FUNCT__
1024a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1037bf11e45SBarry Smith /*@
104564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1057bf11e45SBarry Smith 
10615091d37SBarry Smith     Collective on TS
10715091d37SBarry Smith 
108fb4a63b6SLois Curfman McInnes     Input Parameters:
10915091d37SBarry Smith +   ts - timestep context
11015091d37SBarry Smith -   update - latest solution vector
1117bf11e45SBarry Smith 
112fb4a63b6SLois Curfman McInnes     Output Parameters:
11315091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11415091d37SBarry Smith -   flag - indicates if current timestep was ok
1157bf11e45SBarry Smith 
11615091d37SBarry Smith     Level: advanced
117fee21e36SBarry Smith 
118564e8f4eSLois Curfman McInnes     Notes:
119564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
120564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
121564e8f4eSLois Curfman McInnes 
122564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
123564e8f4eSLois Curfman McInnes 
124564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1257bf11e45SBarry Smith @*/
126ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *flag)
1277bf11e45SBarry Smith {
1287bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
129dfbe8321SBarry Smith   PetscErrorCode ierr;
1307bf11e45SBarry Smith 
1313a40ed3dSBarry Smith   PetscFunctionBegin;
132a7cc72afSBarry Smith   if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);}
1337bf11e45SBarry Smith 
1347bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1357bf11e45SBarry Smith 
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
1377bf11e45SBarry Smith }
1387bf11e45SBarry Smith 
1397bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1407bf11e45SBarry Smith 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
143a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime)
1442d3f70b5SBarry Smith {
1452d3f70b5SBarry Smith   Vec            sol = ts->vec_sol;
146dfbe8321SBarry Smith   PetscErrorCode ierr;
147a7cc72afSBarry Smith   PetscInt       i,max_steps = ts->max_steps,its,lits;
148ace3abfcSBarry Smith   PetscBool      ok;
1497bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
15087828ca2SBarry Smith   PetscReal      current_time_step;
1512d3f70b5SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
1532d3f70b5SBarry Smith   *steps = -ts->steps;
1542d3f70b5SBarry Smith 
1557bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
1562d3f70b5SBarry Smith   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
1577bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
158e6e5fe25SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1597bf11e45SBarry Smith     current_time_step = ts->time_step;
1603f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
161e6e5fe25SBarry Smith     while (PETSC_TRUE) {
1627bf11e45SBarry Smith       ts->ptime  += current_time_step;
163f69a0ea3SMatthew Knepley       ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr);
164b850b91aSLisandro Dalcin       ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1659f954729SBarry Smith       ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
166c9780b6fSBarry Smith       ts->nonlinear_its += its; ts->linear_its += lits;
1677bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
1687bf11e45SBarry Smith       if (ok) break;
1697bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1707bf11e45SBarry Smith       current_time_step = ts->time_step;
1717bf11e45SBarry Smith     }
1727bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
1732d3f70b5SBarry Smith     ts->steps++;
1743f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
1752d3f70b5SBarry Smith   }
176e6e5fe25SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
177e6e5fe25SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
178e6e5fe25SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1792d3f70b5SBarry Smith 
1802d3f70b5SBarry Smith   *steps += ts->steps;
181142b95e3SSatish Balay   *ptime  = ts->ptime;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1832d3f70b5SBarry Smith }
1842d3f70b5SBarry Smith 
1852d3f70b5SBarry Smith /*------------------------------------------------------------*/
1864a2ae208SSatish Balay #undef __FUNCT__
1874a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo"
1886849ba73SBarry Smith static PetscErrorCode TSDestroy_Pseudo(TS ts)
1892d3f70b5SBarry Smith {
1907bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
191dfbe8321SBarry Smith   PetscErrorCode ierr;
1922d3f70b5SBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
194e4ef1970SSatish Balay   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
1957bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
1966f2d6a7bSJed Brown   if (pseudo->xdot) {ierr = VecDestroy(pseudo->xdot);CHKERRQ(ierr);}
197606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1992d3f70b5SBarry Smith }
2002d3f70b5SBarry Smith 
2012d3f70b5SBarry Smith 
2022d3f70b5SBarry Smith /*------------------------------------------------------------*/
2032d3f70b5SBarry Smith 
2044a2ae208SSatish Balay #undef __FUNCT__
2056f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot"
2066f2d6a7bSJed Brown /*
2076f2d6a7bSJed Brown     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
2086f2d6a7bSJed Brown */
2096f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
2102d3f70b5SBarry Smith {
2116f2d6a7bSJed Brown   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
2126f2d6a7bSJed Brown   PetscScalar    mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot;
213dfbe8321SBarry Smith   PetscErrorCode ierr;
214a7cc72afSBarry Smith   PetscInt       i,n;
2152d3f70b5SBarry Smith 
2163a40ed3dSBarry Smith   PetscFunctionBegin;
2176f2d6a7bSJed Brown   ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2186f2d6a7bSJed Brown   ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr);
2196f2d6a7bSJed Brown   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2206f2d6a7bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
2212d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2226f2d6a7bSJed Brown     xdot[i] = mdt*(xnp1[i] - xn[i]);
2232d3f70b5SBarry Smith   }
2246f2d6a7bSJed Brown   ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr);
2256f2d6a7bSJed Brown   ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr);
2266f2d6a7bSJed Brown   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
2276f2d6a7bSJed Brown   *Xdot = pseudo->xdot;
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
2292d3f70b5SBarry Smith }
2302d3f70b5SBarry Smith 
2314a2ae208SSatish Balay #undef __FUNCT__
2320f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo"
2336f2d6a7bSJed Brown /*
2346f2d6a7bSJed Brown     The transient residual is
2356f2d6a7bSJed Brown 
2366f2d6a7bSJed Brown         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
2376f2d6a7bSJed Brown 
2386f2d6a7bSJed Brown     or for ODE,
2396f2d6a7bSJed Brown 
2406f2d6a7bSJed Brown         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
2416f2d6a7bSJed Brown 
2426f2d6a7bSJed Brown     This is the function that must be evaluated for transient simulation and for
2436f2d6a7bSJed Brown     finite difference Jacobians.  On the first Newton step, this algorithm uses
2446f2d6a7bSJed Brown     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
2456f2d6a7bSJed Brown     residual is actually the steady state residual.  Pseudotransient
2466f2d6a7bSJed Brown     continuation as described in the literature is a linearly implicit
2476f2d6a7bSJed Brown     algorithm, it only takes this one Newton step with the steady state
2486f2d6a7bSJed Brown     residual, and then advances to the next time step.
2496f2d6a7bSJed Brown */
2500f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
2512d3f70b5SBarry Smith {
2526f2d6a7bSJed Brown   Vec            Xdot;
253dfbe8321SBarry Smith   PetscErrorCode ierr;
2542d3f70b5SBarry Smith 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2566f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
2576f2d6a7bSJed Brown   ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr);
2586f2d6a7bSJed Brown   PetscFunctionReturn(0);
2596f2d6a7bSJed Brown }
2602d3f70b5SBarry Smith 
2616f2d6a7bSJed Brown #undef __FUNCT__
2620f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
2636f2d6a7bSJed Brown /*
2646f2d6a7bSJed Brown    This constructs the Jacobian needed for SNES.  For DAE, this is
2656f2d6a7bSJed Brown 
2666f2d6a7bSJed Brown        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
2676f2d6a7bSJed Brown 
2686f2d6a7bSJed Brown     and for ODE:
2696f2d6a7bSJed Brown 
2706f2d6a7bSJed Brown        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
2716f2d6a7bSJed Brown */
2720f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
2736f2d6a7bSJed Brown {
2746f2d6a7bSJed Brown   Vec            Xdot;
2756f2d6a7bSJed Brown   PetscErrorCode ierr;
2766f2d6a7bSJed Brown 
2776f2d6a7bSJed Brown   PetscFunctionBegin;
2786f2d6a7bSJed Brown   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
2796f2d6a7bSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr);
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2812d3f70b5SBarry Smith }
2822d3f70b5SBarry Smith 
2832d3f70b5SBarry Smith 
2844a2ae208SSatish Balay #undef __FUNCT__
2854a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2866849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2872d3f70b5SBarry Smith {
2887bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
289dfbe8321SBarry Smith   PetscErrorCode ierr;
2902d3f70b5SBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
2927bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2937bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2946f2d6a7bSJed Brown   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
2950f5c6efeSJed Brown   ierr = SNESSetFunction(ts->snes,pseudo->func,SNESTSFormFunction,ts);CHKERRQ(ierr);
29644c54cacSJed Brown   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
29744c54cacSJed Brown   replace A and we don't want to mess with that.  With -snes_mf, A and B will be replaced as well as the function and
29844c54cacSJed Brown   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
29944c54cacSJed Brown   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
30044c54cacSJed Brown   snes->jacP should be the TS. */
30144c54cacSJed Brown   {
30244c54cacSJed Brown     Mat A,B;
30344c54cacSJed Brown     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
30444c54cacSJed Brown     void *ctx;
30544c54cacSJed Brown     ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr);
30644c54cacSJed Brown     ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr);
30744c54cacSJed Brown   }
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3092d3f70b5SBarry Smith }
3102d3f70b5SBarry Smith /*------------------------------------------------------------*/
3112d3f70b5SBarry Smith 
3124a2ae208SSatish Balay #undef __FUNCT__
313a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
314a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
3152d3f70b5SBarry Smith {
3167bf11e45SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
317dfbe8321SBarry Smith   PetscErrorCode          ierr;
318a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
3192d3f70b5SBarry Smith 
3203a40ed3dSBarry Smith   PetscFunctionBegin;
321a34d58ebSBarry Smith   if (!ctx) {
3227adad957SLisandro Dalcin     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
323a34d58ebSBarry Smith   }
324a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
325a34d58ebSBarry Smith   if (!ctx) {
326a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
327a34d58ebSBarry Smith   }
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
3292d3f70b5SBarry Smith }
3302d3f70b5SBarry Smith 
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
3336849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
3342d3f70b5SBarry Smith {
3354bbc92c1SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
336dfbe8321SBarry Smith   PetscErrorCode          ierr;
337ace3abfcSBarry Smith   PetscBool               flg = PETSC_FALSE;
338a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer;
3392d3f70b5SBarry Smith 
3403a40ed3dSBarry Smith   PetscFunctionBegin;
341b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
342*acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3432d3f70b5SBarry Smith     if (flg) {
3447adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
345a34d58ebSBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
34628aa8177SBarry Smith     }
34790d69ab7SBarry Smith     flg  = PETSC_FALSE;
348*acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
349ca4b7087SBarry Smith     if (flg) {
350ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
351ca4b7087SBarry Smith     }
352419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
353b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3543a40ed3dSBarry Smith   PetscFunctionReturn(0);
3552d3f70b5SBarry Smith }
3562d3f70b5SBarry Smith 
3574a2ae208SSatish Balay #undef __FUNCT__
3584a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3596849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3602d3f70b5SBarry Smith {
3613a40ed3dSBarry Smith   PetscFunctionBegin;
3623a40ed3dSBarry Smith   PetscFunctionReturn(0);
3632d3f70b5SBarry Smith }
3642d3f70b5SBarry Smith 
36582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3664a2ae208SSatish Balay #undef __FUNCT__
3674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
368ac226902SBarry Smith /*@C
36982bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
37082bf6240SBarry Smith    last timestep.
37182bf6240SBarry Smith 
3723f9fe445SBarry Smith    Logically Collective on TS
37315091d37SBarry Smith 
37482bf6240SBarry Smith    Input Parameters:
37515091d37SBarry Smith +  ts - timestep context
37682bf6240SBarry Smith .  dt - user-defined function to verify timestep
37715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
37882bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
37982bf6240SBarry Smith 
38015091d37SBarry Smith    Level: advanced
381fee21e36SBarry Smith 
38282bf6240SBarry Smith    Calling sequence of func:
383ace3abfcSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
38482bf6240SBarry Smith 
38582bf6240SBarry Smith .  update - latest solution vector
38682bf6240SBarry Smith .  ctx - [optional] timestep context
38782bf6240SBarry Smith .  newdt - the timestep to use for the next step
38882bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
38982bf6240SBarry Smith 
39082bf6240SBarry Smith    Notes:
39182bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
39282bf6240SBarry Smith    during the timestepping process.
39382bf6240SBarry Smith 
39482bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
39582bf6240SBarry Smith 
39682bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
39782bf6240SBarry Smith @*/
398ace3abfcSBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
39982bf6240SBarry Smith {
4004ac538c5SBarry Smith   PetscErrorCode ierr;
40182bf6240SBarry Smith 
40282bf6240SBarry Smith   PetscFunctionBegin;
4030700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4044ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
40582bf6240SBarry Smith   PetscFunctionReturn(0);
40682bf6240SBarry Smith }
40782bf6240SBarry Smith 
4084a2ae208SSatish Balay #undef __FUNCT__
4094a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
41082bf6240SBarry Smith /*@
41182bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
41282bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
41382bf6240SBarry Smith 
4143f9fe445SBarry Smith    Logically Collective on TS
415fee21e36SBarry Smith 
41615091d37SBarry Smith     Input Parameters:
41715091d37SBarry Smith +   ts - the timestep context
41815091d37SBarry Smith -   inc - the scaling factor >= 1.0
41915091d37SBarry Smith 
42082bf6240SBarry Smith     Options Database Key:
42182bf6240SBarry Smith $    -ts_pseudo_increment <increment>
42282bf6240SBarry Smith 
42315091d37SBarry Smith     Level: advanced
42415091d37SBarry Smith 
42582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42682bf6240SBarry Smith 
42782bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42882bf6240SBarry Smith @*/
42963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
43082bf6240SBarry Smith {
4314ac538c5SBarry Smith   PetscErrorCode ierr;
43282bf6240SBarry Smith 
43382bf6240SBarry Smith   PetscFunctionBegin;
4340700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
435c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,inc,2);
4364ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
43782bf6240SBarry Smith   PetscFunctionReturn(0);
43882bf6240SBarry Smith }
43982bf6240SBarry Smith 
4404a2ae208SSatish Balay #undef __FUNCT__
4414a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
44282bf6240SBarry Smith /*@
44382bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
44482bf6240SBarry Smith     is computed via the formula
44582bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
44682bf6240SBarry Smith       rather than the default update,
44782bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
44882bf6240SBarry Smith 
4493f9fe445SBarry Smith    Logically Collective on TS
45015091d37SBarry Smith 
45182bf6240SBarry Smith     Input Parameter:
45282bf6240SBarry Smith .   ts - the timestep context
45382bf6240SBarry Smith 
45482bf6240SBarry Smith     Options Database Key:
45582bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
45682bf6240SBarry Smith 
45715091d37SBarry Smith     Level: advanced
45815091d37SBarry Smith 
45982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
46082bf6240SBarry Smith 
46182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
46282bf6240SBarry Smith @*/
46363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt(TS ts)
46482bf6240SBarry Smith {
4654ac538c5SBarry Smith   PetscErrorCode ierr;
46682bf6240SBarry Smith 
46782bf6240SBarry Smith   PetscFunctionBegin;
4680700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4694ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
47082bf6240SBarry Smith   PetscFunctionReturn(0);
47182bf6240SBarry Smith }
47282bf6240SBarry Smith 
47382bf6240SBarry Smith 
4744a2ae208SSatish Balay #undef __FUNCT__
4754a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
476ac226902SBarry Smith /*@C
47782bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
47882bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
47982bf6240SBarry Smith 
4803f9fe445SBarry Smith    Logically Collective on TS
48115091d37SBarry Smith 
48282bf6240SBarry Smith    Input Parameters:
48315091d37SBarry Smith +  ts - timestep context
48482bf6240SBarry Smith .  dt - function to compute timestep
48515091d37SBarry Smith -  ctx - [optional] user-defined context for private data
48682bf6240SBarry Smith          required by the function (may be PETSC_NULL)
48782bf6240SBarry Smith 
48815091d37SBarry Smith    Level: intermediate
489fee21e36SBarry Smith 
49082bf6240SBarry Smith    Calling sequence of func:
49187828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
49282bf6240SBarry Smith 
49382bf6240SBarry Smith .  newdt - the newly computed timestep
49482bf6240SBarry Smith .  ctx - [optional] timestep context
49582bf6240SBarry Smith 
49682bf6240SBarry Smith    Notes:
49782bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
49882bf6240SBarry Smith    during the timestepping process.
49982bf6240SBarry Smith 
50082bf6240SBarry Smith .keywords: timestep, pseudo, set
50182bf6240SBarry Smith 
50282bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
50382bf6240SBarry Smith @*/
50463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
50582bf6240SBarry Smith {
5064ac538c5SBarry Smith   PetscErrorCode ierr;
50782bf6240SBarry Smith 
50882bf6240SBarry Smith   PetscFunctionBegin;
5090700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5104ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
51182bf6240SBarry Smith   PetscFunctionReturn(0);
51282bf6240SBarry Smith }
51382bf6240SBarry Smith 
51482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
51582bf6240SBarry Smith 
516ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
517fb2e594dSBarry Smith EXTERN_C_BEGIN
5184a2ae208SSatish Balay #undef __FUNCT__
5194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
52063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
52182bf6240SBarry Smith {
52282bf6240SBarry Smith   TS_Pseudo *pseudo;
52382bf6240SBarry Smith 
52482bf6240SBarry Smith   PetscFunctionBegin;
52582bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
52682bf6240SBarry Smith   pseudo->verify      = dt;
52782bf6240SBarry Smith   pseudo->verifyctx   = ctx;
52882bf6240SBarry Smith   PetscFunctionReturn(0);
52982bf6240SBarry Smith }
530fb2e594dSBarry Smith EXTERN_C_END
53182bf6240SBarry Smith 
532fb2e594dSBarry Smith EXTERN_C_BEGIN
5334a2ae208SSatish Balay #undef __FUNCT__
5344a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
53563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
53682bf6240SBarry Smith {
5374bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53882bf6240SBarry Smith 
53982bf6240SBarry Smith   PetscFunctionBegin;
54082bf6240SBarry Smith   pseudo->dt_increment = inc;
54182bf6240SBarry Smith   PetscFunctionReturn(0);
54282bf6240SBarry Smith }
543fb2e594dSBarry Smith EXTERN_C_END
54482bf6240SBarry Smith 
545fb2e594dSBarry Smith EXTERN_C_BEGIN
5464a2ae208SSatish Balay #undef __FUNCT__
5474a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
54863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
54982bf6240SBarry Smith {
5504bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55182bf6240SBarry Smith 
55282bf6240SBarry Smith   PetscFunctionBegin;
5534bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
55482bf6240SBarry Smith   PetscFunctionReturn(0);
55582bf6240SBarry Smith }
556fb2e594dSBarry Smith EXTERN_C_END
55782bf6240SBarry Smith 
5586849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
559fb2e594dSBarry Smith EXTERN_C_BEGIN
5604a2ae208SSatish Balay #undef __FUNCT__
5614a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
56263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
56382bf6240SBarry Smith {
5644bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
56582bf6240SBarry Smith 
56682bf6240SBarry Smith   PetscFunctionBegin;
56782bf6240SBarry Smith   pseudo->dt      = dt;
56882bf6240SBarry Smith   pseudo->dtctx   = ctx;
56982bf6240SBarry Smith   PetscFunctionReturn(0);
57082bf6240SBarry Smith }
571fb2e594dSBarry Smith EXTERN_C_END
57282bf6240SBarry Smith 
57382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
57410e6a065SJed Brown /*MC
57510e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
57682bf6240SBarry Smith 
57710e6a065SJed Brown   This method solves equations of the form
57810e6a065SJed Brown 
57910e6a065SJed Brown $    F(X,Xdot) = 0
58010e6a065SJed Brown 
58110e6a065SJed Brown   for steady state using the iteration
58210e6a065SJed Brown 
58310e6a065SJed Brown $    [G'] S = -F(X,0)
58410e6a065SJed Brown $    X += S
58510e6a065SJed Brown 
58610e6a065SJed Brown   where
58710e6a065SJed Brown 
58810e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
58910e6a065SJed Brown 
5906f2d6a7bSJed Brown   This is linearly-implicit Euler with the residual always evaluated "at steady
5916f2d6a7bSJed Brown   state".  See note below.
59210e6a065SJed Brown 
59310e6a065SJed Brown   Options database keys:
59410e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
59510e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
59610e6a065SJed Brown 
59710e6a065SJed Brown   Level: beginner
59810e6a065SJed Brown 
59910e6a065SJed Brown   References:
60010e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
60110e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
60210e6a065SJed Brown 
60310e6a065SJed Brown   Notes:
6046f2d6a7bSJed Brown   The residual computed by this method includes the transient term (Xdot is computed instead of
6056f2d6a7bSJed Brown   always being zero), but since the prediction from the last step is always the solution from the
6066f2d6a7bSJed Brown   last step, on the first Newton iteration we have
6076f2d6a7bSJed Brown 
6086f2d6a7bSJed Brown $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
6096f2d6a7bSJed Brown 
6106f2d6a7bSJed Brown   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
6116f2d6a7bSJed Brown   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
6126f2d6a7bSJed Brown   algorithm is no longer the one described in the referenced papers.
61310e6a065SJed Brown 
61410e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
61510e6a065SJed Brown 
61610e6a065SJed Brown M*/
617fb2e594dSBarry Smith EXTERN_C_BEGIN
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
62063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Pseudo(TS ts)
6212d3f70b5SBarry Smith {
6227bf11e45SBarry Smith   TS_Pseudo      *pseudo;
623dfbe8321SBarry Smith   PetscErrorCode ierr;
6242d3f70b5SBarry Smith 
6253a40ed3dSBarry Smith   PetscFunctionBegin;
626000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
627000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
6282d3f70b5SBarry Smith 
62917186662SBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
630000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
631000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
632000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6330f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
6340f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
6357bf11e45SBarry Smith 
6367bf11e45SBarry Smith   /* create the required nonlinear solver context */
6377adad957SLisandro Dalcin   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
63817efb626SBarry Smith   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
6392d3f70b5SBarry Smith 
64038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6417bf11e45SBarry Smith   ts->data = (void*)pseudo;
6422d3f70b5SBarry Smith 
64328aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6444bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
64528aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
64682bf6240SBarry Smith 
647f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
648e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6490c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
650f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
651e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6520c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
653f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
654e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6550c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6560c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6570c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6580c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6593a40ed3dSBarry Smith   PetscFunctionReturn(0);
6602d3f70b5SBarry Smith }
661fb2e594dSBarry Smith EXTERN_C_END
6622d3f70b5SBarry Smith 
6634a2ae208SSatish Balay #undef __FUNCT__
6644a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
66582bf6240SBarry Smith /*@C
66682bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
66782bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
66828aa8177SBarry Smith 
66915091d37SBarry Smith    Collective on TS
67015091d37SBarry Smith 
67128aa8177SBarry Smith    Input Parameters:
67228aa8177SBarry Smith .  ts - the timestep context
67382bf6240SBarry Smith .  dtctx - unused timestep context
67428aa8177SBarry Smith 
67582bf6240SBarry Smith    Output Parameter:
67682bf6240SBarry Smith .  newdt - the timestep to use for the next step
67728aa8177SBarry Smith 
67815091d37SBarry Smith    Level: advanced
67915091d37SBarry Smith 
68082bf6240SBarry Smith .keywords: timestep, pseudo, default
681564e8f4eSLois Curfman McInnes 
68282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
68328aa8177SBarry Smith @*/
68463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
68528aa8177SBarry Smith {
68682bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
68787828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
688dfbe8321SBarry Smith   PetscErrorCode ierr;
68928aa8177SBarry Smith 
6903a40ed3dSBarry Smith   PetscFunctionBegin;
69182bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
69282bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
69382bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
69482bf6240SBarry Smith     /* first time through so compute initial function norm */
69582bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
69682bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
69782bf6240SBarry Smith   }
69882bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
69982bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
700004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
70182bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
70282bf6240SBarry Smith   } else {
70382bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
70482bf6240SBarry Smith   }
70582bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
7063a40ed3dSBarry Smith   PetscFunctionReturn(0);
70728aa8177SBarry Smith }
708