xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 10e6a065829bb40709f4177632bce2eca368f9ba)
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 */
112d3f70b5SBarry Smith   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
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;
17a7cc72afSBarry Smith   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscTruth*); /* 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 */
244bbc92c1SBarry Smith   PetscTruth 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 @*/
9363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *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 @*/
12663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *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;
148a7cc72afSBarry Smith   PetscTruth     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);}
1967bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
197606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1992d3f70b5SBarry Smith }
2002d3f70b5SBarry Smith 
2012d3f70b5SBarry Smith 
2022d3f70b5SBarry Smith /*------------------------------------------------------------*/
2032d3f70b5SBarry Smith 
2042d3f70b5SBarry Smith /*
2052d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
2062d3f70b5SBarry Smith 
2072d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
2082d3f70b5SBarry Smith */
2094a2ae208SSatish Balay #undef __FUNCT__
2104a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction"
211dfbe8321SBarry Smith PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2122d3f70b5SBarry Smith {
2132d3f70b5SBarry Smith   TS             ts = (TS) ctx;
214ea709b57SSatish Balay   PetscScalar    mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
215dfbe8321SBarry Smith   PetscErrorCode ierr;
216a7cc72afSBarry Smith   PetscInt       i,n;
2172d3f70b5SBarry Smith 
2183a40ed3dSBarry Smith   PetscFunctionBegin;
2192d3f70b5SBarry Smith   /* apply user provided function */
2202d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
2217bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2222d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
2232d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
2242d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
2252d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
2262d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2272d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2282d3f70b5SBarry Smith   }
229888f2ed8SSatish Balay   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
230888f2ed8SSatish Balay   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
231888f2ed8SSatish Balay   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
2323a40ed3dSBarry Smith   PetscFunctionReturn(0);
2332d3f70b5SBarry Smith }
2342d3f70b5SBarry Smith 
2352d3f70b5SBarry Smith /*
2362d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2372d3f70b5SBarry Smith 
2382d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2392d3f70b5SBarry Smith */
2404a2ae208SSatish Balay #undef __FUNCT__
2414a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian"
242dfbe8321SBarry Smith PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2432d3f70b5SBarry Smith {
2442d3f70b5SBarry Smith   TS             ts = (TS) ctx;
245dfbe8321SBarry Smith   PetscErrorCode ierr;
2462d3f70b5SBarry Smith 
2473a40ed3dSBarry Smith   PetscFunctionBegin;
2482d3f70b5SBarry Smith   /* construct users Jacobian */
249a7a1495cSBarry Smith   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
2502d3f70b5SBarry Smith 
251614f654bSBarry Smith   /* shift and scale Jacobian */
252b0f6e734SBarry Smith   ierr = TSScaleShiftMatrices(ts,*AA,*BB,*str);CHKERRQ(ierr);
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2542d3f70b5SBarry Smith }
2552d3f70b5SBarry Smith 
2562d3f70b5SBarry Smith 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2596849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2602d3f70b5SBarry Smith {
2617bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
262dfbe8321SBarry Smith   PetscErrorCode ierr;
2632d3f70b5SBarry Smith 
2643a40ed3dSBarry Smith   PetscFunctionBegin;
2657bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2667bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2677bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2688beb423aSHong Zhang   ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
2693a40ed3dSBarry Smith   PetscFunctionReturn(0);
2702d3f70b5SBarry Smith }
2712d3f70b5SBarry Smith /*------------------------------------------------------------*/
2722d3f70b5SBarry Smith 
2734a2ae208SSatish Balay #undef __FUNCT__
274a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault"
275a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
2762d3f70b5SBarry Smith {
2777bf11e45SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
278dfbe8321SBarry Smith   PetscErrorCode          ierr;
279a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
2802d3f70b5SBarry Smith 
2813a40ed3dSBarry Smith   PetscFunctionBegin;
282a34d58ebSBarry Smith   if (!ctx) {
2837adad957SLisandro Dalcin     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
284a34d58ebSBarry Smith   }
285a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
286a34d58ebSBarry Smith   if (!ctx) {
287a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
288a34d58ebSBarry Smith   }
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
2902d3f70b5SBarry Smith }
2912d3f70b5SBarry Smith 
2924a2ae208SSatish Balay #undef __FUNCT__
2934a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
2946849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
2952d3f70b5SBarry Smith {
2964bbc92c1SBarry Smith   TS_Pseudo               *pseudo = (TS_Pseudo*)ts->data;
297dfbe8321SBarry Smith   PetscErrorCode          ierr;
29890d69ab7SBarry Smith   PetscTruth              flg = PETSC_FALSE;
299a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer;
3002d3f70b5SBarry Smith 
3013a40ed3dSBarry Smith   PetscFunctionBegin;
302b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
30390d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3042d3f70b5SBarry Smith     if (flg) {
3057adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
306a34d58ebSBarry Smith       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
30728aa8177SBarry Smith     }
30890d69ab7SBarry Smith     flg  = PETSC_FALSE;
30990d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
310ca4b7087SBarry Smith     if (flg) {
311ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
312ca4b7087SBarry Smith     }
313419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
314b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
3162d3f70b5SBarry Smith }
3172d3f70b5SBarry Smith 
3184a2ae208SSatish Balay #undef __FUNCT__
3194a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3206849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3212d3f70b5SBarry Smith {
3223a40ed3dSBarry Smith   PetscFunctionBegin;
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
3242d3f70b5SBarry Smith }
3252d3f70b5SBarry Smith 
32682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3274a2ae208SSatish Balay #undef __FUNCT__
3284a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
329ac226902SBarry Smith /*@C
33082bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
33182bf6240SBarry Smith    last timestep.
33282bf6240SBarry Smith 
33315091d37SBarry Smith    Collective on TS
33415091d37SBarry Smith 
33582bf6240SBarry Smith    Input Parameters:
33615091d37SBarry Smith +  ts - timestep context
33782bf6240SBarry Smith .  dt - user-defined function to verify timestep
33815091d37SBarry Smith -  ctx - [optional] user-defined context for private data
33982bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
34082bf6240SBarry Smith 
34115091d37SBarry Smith    Level: advanced
342fee21e36SBarry Smith 
34382bf6240SBarry Smith    Calling sequence of func:
344a7cc72afSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag);
34582bf6240SBarry Smith 
34682bf6240SBarry Smith .  update - latest solution vector
34782bf6240SBarry Smith .  ctx - [optional] timestep context
34882bf6240SBarry Smith .  newdt - the timestep to use for the next step
34982bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
35082bf6240SBarry Smith 
35182bf6240SBarry Smith    Notes:
35282bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
35382bf6240SBarry Smith    during the timestepping process.
35482bf6240SBarry Smith 
35582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
35682bf6240SBarry Smith 
35782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
35882bf6240SBarry Smith @*/
35963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx)
36082bf6240SBarry Smith {
361a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *);
36282bf6240SBarry Smith 
36382bf6240SBarry Smith   PetscFunctionBegin;
3644482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
36582bf6240SBarry Smith 
366c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
36782bf6240SBarry Smith   if (f) {
36882bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
36982bf6240SBarry Smith   }
37082bf6240SBarry Smith   PetscFunctionReturn(0);
37182bf6240SBarry Smith }
37282bf6240SBarry Smith 
3734a2ae208SSatish Balay #undef __FUNCT__
3744a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
37582bf6240SBarry Smith /*@
37682bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
37782bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
37882bf6240SBarry Smith 
379fee21e36SBarry Smith    Collective on TS
380fee21e36SBarry Smith 
38115091d37SBarry Smith     Input Parameters:
38215091d37SBarry Smith +   ts - the timestep context
38315091d37SBarry Smith -   inc - the scaling factor >= 1.0
38415091d37SBarry Smith 
38582bf6240SBarry Smith     Options Database Key:
38682bf6240SBarry Smith $    -ts_pseudo_increment <increment>
38782bf6240SBarry Smith 
38815091d37SBarry Smith     Level: advanced
38915091d37SBarry Smith 
39082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
39182bf6240SBarry Smith 
39282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
39382bf6240SBarry Smith @*/
39463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
39582bf6240SBarry Smith {
396dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscReal);
39782bf6240SBarry Smith 
39882bf6240SBarry Smith   PetscFunctionBegin;
3994482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
40082bf6240SBarry Smith 
401c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
40282bf6240SBarry Smith   if (f) {
40382bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
40482bf6240SBarry Smith   }
40582bf6240SBarry Smith   PetscFunctionReturn(0);
40682bf6240SBarry Smith }
40782bf6240SBarry Smith 
4084a2ae208SSatish Balay #undef __FUNCT__
4094a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
41082bf6240SBarry Smith /*@
41182bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
41282bf6240SBarry Smith     is computed via the formula
41382bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
41482bf6240SBarry Smith       rather than the default update,
41582bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
41682bf6240SBarry Smith 
41715091d37SBarry Smith    Collective on TS
41815091d37SBarry Smith 
41982bf6240SBarry Smith     Input Parameter:
42082bf6240SBarry Smith .   ts - the timestep context
42182bf6240SBarry Smith 
42282bf6240SBarry Smith     Options Database Key:
42382bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
42482bf6240SBarry Smith 
42515091d37SBarry Smith     Level: advanced
42615091d37SBarry Smith 
42782bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42882bf6240SBarry Smith 
42982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
43082bf6240SBarry Smith @*/
43163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt(TS ts)
43282bf6240SBarry Smith {
433dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS);
43482bf6240SBarry Smith 
43582bf6240SBarry Smith   PetscFunctionBegin;
4364482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
43782bf6240SBarry Smith 
438c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
43982bf6240SBarry Smith   if (f) {
44082bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
44182bf6240SBarry Smith   }
44282bf6240SBarry Smith   PetscFunctionReturn(0);
44382bf6240SBarry Smith }
44482bf6240SBarry Smith 
44582bf6240SBarry Smith 
4464a2ae208SSatish Balay #undef __FUNCT__
4474a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
448ac226902SBarry Smith /*@C
44982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
45082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
45182bf6240SBarry Smith 
45215091d37SBarry Smith    Collective on TS
45315091d37SBarry Smith 
45482bf6240SBarry Smith    Input Parameters:
45515091d37SBarry Smith +  ts - timestep context
45682bf6240SBarry Smith .  dt - function to compute timestep
45715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
45882bf6240SBarry Smith          required by the function (may be PETSC_NULL)
45982bf6240SBarry Smith 
46015091d37SBarry Smith    Level: intermediate
461fee21e36SBarry Smith 
46282bf6240SBarry Smith    Calling sequence of func:
46387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
46482bf6240SBarry Smith 
46582bf6240SBarry Smith .  newdt - the newly computed timestep
46682bf6240SBarry Smith .  ctx - [optional] timestep context
46782bf6240SBarry Smith 
46882bf6240SBarry Smith    Notes:
46982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
47082bf6240SBarry Smith    during the timestepping process.
47182bf6240SBarry Smith 
47282bf6240SBarry Smith .keywords: timestep, pseudo, set
47382bf6240SBarry Smith 
47482bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
47582bf6240SBarry Smith @*/
47663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
47782bf6240SBarry Smith {
4786849ba73SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *);
47982bf6240SBarry Smith 
48082bf6240SBarry Smith   PetscFunctionBegin;
4814482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
48282bf6240SBarry Smith 
483c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
48482bf6240SBarry Smith   if (f) {
48582bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
48682bf6240SBarry Smith   }
48782bf6240SBarry Smith   PetscFunctionReturn(0);
48882bf6240SBarry Smith }
48982bf6240SBarry Smith 
49082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
49182bf6240SBarry Smith 
492a7cc72afSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
493fb2e594dSBarry Smith EXTERN_C_BEGIN
4944a2ae208SSatish Balay #undef __FUNCT__
4954a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
49663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
49782bf6240SBarry Smith {
49882bf6240SBarry Smith   TS_Pseudo *pseudo;
49982bf6240SBarry Smith 
50082bf6240SBarry Smith   PetscFunctionBegin;
50182bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
50282bf6240SBarry Smith   pseudo->verify      = dt;
50382bf6240SBarry Smith   pseudo->verifyctx   = ctx;
50482bf6240SBarry Smith   PetscFunctionReturn(0);
50582bf6240SBarry Smith }
506fb2e594dSBarry Smith EXTERN_C_END
50782bf6240SBarry Smith 
508fb2e594dSBarry Smith EXTERN_C_BEGIN
5094a2ae208SSatish Balay #undef __FUNCT__
5104a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
51163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
51282bf6240SBarry Smith {
5134bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
51482bf6240SBarry Smith 
51582bf6240SBarry Smith   PetscFunctionBegin;
51682bf6240SBarry Smith   pseudo->dt_increment = inc;
51782bf6240SBarry Smith   PetscFunctionReturn(0);
51882bf6240SBarry Smith }
519fb2e594dSBarry Smith EXTERN_C_END
52082bf6240SBarry Smith 
521fb2e594dSBarry Smith EXTERN_C_BEGIN
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
52463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
52582bf6240SBarry Smith {
5264bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
52782bf6240SBarry Smith 
52882bf6240SBarry Smith   PetscFunctionBegin;
5294bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
53082bf6240SBarry Smith   PetscFunctionReturn(0);
53182bf6240SBarry Smith }
532fb2e594dSBarry Smith EXTERN_C_END
53382bf6240SBarry Smith 
5346849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
535fb2e594dSBarry Smith EXTERN_C_BEGIN
5364a2ae208SSatish Balay #undef __FUNCT__
5374a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
53863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
53982bf6240SBarry Smith {
5404bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54182bf6240SBarry Smith 
54282bf6240SBarry Smith   PetscFunctionBegin;
54382bf6240SBarry Smith   pseudo->dt      = dt;
54482bf6240SBarry Smith   pseudo->dtctx   = ctx;
54582bf6240SBarry Smith   PetscFunctionReturn(0);
54682bf6240SBarry Smith }
547fb2e594dSBarry Smith EXTERN_C_END
54882bf6240SBarry Smith 
54982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
550*10e6a065SJed Brown /*MC
551*10e6a065SJed Brown       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
55282bf6240SBarry Smith 
553*10e6a065SJed Brown   This method solves equations of the form
554*10e6a065SJed Brown 
555*10e6a065SJed Brown $    F(X,Xdot) = 0
556*10e6a065SJed Brown 
557*10e6a065SJed Brown   for steady state using the iteration
558*10e6a065SJed Brown 
559*10e6a065SJed Brown $    [G'] S = -F(X,0)
560*10e6a065SJed Brown $    X += S
561*10e6a065SJed Brown 
562*10e6a065SJed Brown   where
563*10e6a065SJed Brown 
564*10e6a065SJed Brown $    G(Y) = F(Y,(Y-X)/dt)
565*10e6a065SJed Brown 
566*10e6a065SJed Brown   This is very much like linearly-implicit Euler, except that the residual is evaluated "at steady
567*10e6a065SJed Brown   state".
568*10e6a065SJed Brown 
569*10e6a065SJed Brown   Options database keys:
570*10e6a065SJed Brown +  -ts_pseudo_increment <real> - ratio of increase dt
571*10e6a065SJed Brown -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
572*10e6a065SJed Brown 
573*10e6a065SJed Brown   Level: beginner
574*10e6a065SJed Brown 
575*10e6a065SJed Brown   References:
576*10e6a065SJed Brown   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
577*10e6a065SJed Brown   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
578*10e6a065SJed Brown 
579*10e6a065SJed Brown   Notes:
580*10e6a065SJed Brown   At the moment, this uses uses a transient residual, so it is not the pseudotransient continuation described in these papers.
581*10e6a065SJed Brown 
582*10e6a065SJed Brown .seealso:  TSCreate(), TS, TSSetType()
583*10e6a065SJed Brown 
584*10e6a065SJed Brown M*/
585fb2e594dSBarry Smith EXTERN_C_BEGIN
5864a2ae208SSatish Balay #undef __FUNCT__
5874a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
58863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Pseudo(TS ts)
5892d3f70b5SBarry Smith {
5907bf11e45SBarry Smith   TS_Pseudo      *pseudo;
591dfbe8321SBarry Smith   PetscErrorCode ierr;
5922d3f70b5SBarry Smith 
5933a40ed3dSBarry Smith   PetscFunctionBegin;
594000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
595000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
5962d3f70b5SBarry Smith 
5972d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
59829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
5992d3f70b5SBarry Smith   }
6008beb423aSHong Zhang   if (!ts->Arhs) {
60129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
6022d3f70b5SBarry Smith   }
603273d9f13SBarry Smith 
604000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
605000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
606000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
6077bf11e45SBarry Smith 
6087bf11e45SBarry Smith   /* create the required nonlinear solver context */
6097adad957SLisandro Dalcin   ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
61017efb626SBarry Smith   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
6112d3f70b5SBarry Smith 
61238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
6137bf11e45SBarry Smith   ts->data = (void*)pseudo;
6142d3f70b5SBarry Smith 
61528aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
6164bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
61728aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
61882bf6240SBarry Smith 
619f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
620e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
6210c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
622f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
623e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
6240c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
625f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
626e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
6270c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
6280c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
6290c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
6300c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6313a40ed3dSBarry Smith   PetscFunctionReturn(0);
6322d3f70b5SBarry Smith }
633fb2e594dSBarry Smith EXTERN_C_END
6342d3f70b5SBarry Smith 
6354a2ae208SSatish Balay #undef __FUNCT__
6364a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
63782bf6240SBarry Smith /*@C
63882bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
63982bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
64028aa8177SBarry Smith 
64115091d37SBarry Smith    Collective on TS
64215091d37SBarry Smith 
64328aa8177SBarry Smith    Input Parameters:
64428aa8177SBarry Smith .  ts - the timestep context
64582bf6240SBarry Smith .  dtctx - unused timestep context
64628aa8177SBarry Smith 
64782bf6240SBarry Smith    Output Parameter:
64882bf6240SBarry Smith .  newdt - the timestep to use for the next step
64928aa8177SBarry Smith 
65015091d37SBarry Smith    Level: advanced
65115091d37SBarry Smith 
65282bf6240SBarry Smith .keywords: timestep, pseudo, default
653564e8f4eSLois Curfman McInnes 
65482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
65528aa8177SBarry Smith @*/
65663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
65728aa8177SBarry Smith {
65882bf6240SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
65987828ca2SBarry Smith   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
660dfbe8321SBarry Smith   PetscErrorCode ierr;
66128aa8177SBarry Smith 
6623a40ed3dSBarry Smith   PetscFunctionBegin;
66382bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
66482bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
66582bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
66682bf6240SBarry Smith     /* first time through so compute initial function norm */
66782bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
66882bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
66982bf6240SBarry Smith   }
67082bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
67182bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
672004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
67382bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
67482bf6240SBarry Smith   } else {
67582bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
67682bf6240SBarry Smith   }
67782bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6783a40ed3dSBarry Smith   PetscFunctionReturn(0);
67928aa8177SBarry Smith }
680