xref: /petsc/src/ts/impls/pseudo/posindep.c (revision b0f6e734dc4c6779e1d44bf1592310611e4e1d43)
12d3f70b5SBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
32d3f70b5SBarry Smith */
4e090d566SSatish Balay #include "src/ts/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 */
92d3f70b5SBarry Smith   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
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;
15a7cc72afSBarry Smith   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscTruth*); /* verify previous timestep and related context */
167bf11e45SBarry Smith   void           *verifyctx;
172d3f70b5SBarry Smith 
1887828ca2SBarry Smith   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
1987828ca2SBarry Smith   PetscReal  fnorm_previous;
2028aa8177SBarry Smith 
2187828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
224bbc92c1SBarry Smith   PetscTruth increment_dt_from_initial_dt;
237bf11e45SBarry Smith } TS_Pseudo;
242d3f70b5SBarry Smith 
252d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
262d3f70b5SBarry Smith 
274a2ae208SSatish Balay #undef __FUNCT__
284a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
297bf11e45SBarry Smith /*@
307bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
31564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
322d3f70b5SBarry Smith 
3315091d37SBarry Smith     Collective on TS
3415091d37SBarry Smith 
357bf11e45SBarry Smith     Input Parameter:
367bf11e45SBarry Smith .   ts - timestep context
377bf11e45SBarry Smith 
387bf11e45SBarry Smith     Output Parameter:
39fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
40fb4a63b6SLois Curfman McInnes 
4115091d37SBarry Smith     Level: advanced
42564e8f4eSLois Curfman McInnes 
43564e8f4eSLois Curfman McInnes     Notes:
44564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
45564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
46564e8f4eSLois Curfman McInnes 
47fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
48564e8f4eSLois Curfman McInnes 
49564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
507bf11e45SBarry Smith @*/
51dfbe8321SBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
527bf11e45SBarry Smith {
537bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
557bf11e45SBarry Smith 
563a40ed3dSBarry Smith   PetscFunctionBegin;
57d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
587bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
59d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
617bf11e45SBarry Smith }
627bf11e45SBarry Smith 
637bf11e45SBarry Smith 
647bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
654a2ae208SSatish Balay #undef __FUNCT__
664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
677bf11e45SBarry Smith /*@C
68639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
697bf11e45SBarry Smith 
7015091d37SBarry Smith    Collective on TS
7115091d37SBarry Smith 
727bf11e45SBarry Smith    Input Parameters:
7315091d37SBarry Smith +  ts - the timestep context
747bf11e45SBarry Smith .  dtctx - unused timestep context
7515091d37SBarry Smith -  update - latest solution vector
767bf11e45SBarry Smith 
77564e8f4eSLois Curfman McInnes    Output Parameters:
7815091d37SBarry Smith +  newdt - the timestep to use for the next step
7915091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
807bf11e45SBarry Smith 
8115091d37SBarry Smith    Level: advanced
82fee21e36SBarry Smith 
83564e8f4eSLois Curfman McInnes    Note:
84564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
85564e8f4eSLois Curfman McInnes    timestep.
86564e8f4eSLois Curfman McInnes 
87564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
88564e8f4eSLois Curfman McInnes 
89564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
907bf11e45SBarry Smith @*/
91a7cc72afSBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscTruth *flag)
927bf11e45SBarry Smith {
933a40ed3dSBarry Smith   PetscFunctionBegin;
94a7cc72afSBarry Smith   *flag = PETSC_TRUE;
953a40ed3dSBarry Smith   PetscFunctionReturn(0);
967bf11e45SBarry Smith }
977bf11e45SBarry Smith 
987bf11e45SBarry Smith 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1017bf11e45SBarry Smith /*@
102564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1037bf11e45SBarry Smith 
10415091d37SBarry Smith     Collective on TS
10515091d37SBarry Smith 
106fb4a63b6SLois Curfman McInnes     Input Parameters:
10715091d37SBarry Smith +   ts - timestep context
10815091d37SBarry Smith -   update - latest solution vector
1097bf11e45SBarry Smith 
110fb4a63b6SLois Curfman McInnes     Output Parameters:
11115091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11215091d37SBarry Smith -   flag - indicates if current timestep was ok
1137bf11e45SBarry Smith 
11415091d37SBarry Smith     Level: advanced
115fee21e36SBarry Smith 
116564e8f4eSLois Curfman McInnes     Notes:
117564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
118564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
119564e8f4eSLois Curfman McInnes 
120564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
121564e8f4eSLois Curfman McInnes 
122564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1237bf11e45SBarry Smith @*/
124a7cc72afSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscTruth *flag)
1257bf11e45SBarry Smith {
1267bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
127dfbe8321SBarry Smith   PetscErrorCode ierr;
1287bf11e45SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionBegin;
130a7cc72afSBarry Smith   if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);}
1317bf11e45SBarry Smith 
1327bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1337bf11e45SBarry Smith 
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1357bf11e45SBarry Smith }
1367bf11e45SBarry Smith 
1377bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1387bf11e45SBarry Smith 
1394a2ae208SSatish Balay #undef __FUNCT__
1404a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
141a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime)
1422d3f70b5SBarry Smith {
1432d3f70b5SBarry Smith   Vec            sol = ts->vec_sol;
144dfbe8321SBarry Smith   PetscErrorCode ierr;
145a7cc72afSBarry Smith   PetscInt       i,max_steps = ts->max_steps,its,lits;
146a7cc72afSBarry Smith   PetscTruth     ok;
1477bf11e45SBarry Smith   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
14887828ca2SBarry Smith   PetscReal      current_time_step;
1492d3f70b5SBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1512d3f70b5SBarry Smith   *steps = -ts->steps;
1522d3f70b5SBarry Smith 
1537bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
1542d3f70b5SBarry Smith   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
1557bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
156e6e5fe25SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1577bf11e45SBarry Smith     current_time_step = ts->time_step;
158e6e5fe25SBarry Smith     while (PETSC_TRUE) {
1597bf11e45SBarry Smith       ts->ptime  += current_time_step;
160c9780b6fSBarry Smith       ierr = SNESSolve(ts->snes,pseudo->update);CHKERRQ(ierr);
161f2267985SLois Curfman McInnes       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
1629f954729SBarry Smith       ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
163c9780b6fSBarry Smith       ts->nonlinear_its += its; ts->linear_its += lits;
1647bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
1657bf11e45SBarry Smith       if (ok) break;
1667bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1677bf11e45SBarry Smith       current_time_step = ts->time_step;
1687bf11e45SBarry Smith     }
1697bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
1702d3f70b5SBarry Smith     ts->steps++;
1712d3f70b5SBarry Smith   }
172e6e5fe25SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
173e6e5fe25SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
174e6e5fe25SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1752d3f70b5SBarry Smith 
1762d3f70b5SBarry Smith   *steps += ts->steps;
177142b95e3SSatish Balay   *ptime  = ts->ptime;
1783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1792d3f70b5SBarry Smith }
1802d3f70b5SBarry Smith 
1812d3f70b5SBarry Smith /*------------------------------------------------------------*/
1824a2ae208SSatish Balay #undef __FUNCT__
1834a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo"
1846849ba73SBarry Smith static PetscErrorCode TSDestroy_Pseudo(TS ts)
1852d3f70b5SBarry Smith {
1867bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
187dfbe8321SBarry Smith   PetscErrorCode ierr;
1882d3f70b5SBarry Smith 
1893a40ed3dSBarry Smith   PetscFunctionBegin;
190e4ef1970SSatish Balay   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
1917bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
1927bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
193606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1943a40ed3dSBarry Smith   PetscFunctionReturn(0);
1952d3f70b5SBarry Smith }
1962d3f70b5SBarry Smith 
1972d3f70b5SBarry Smith 
1982d3f70b5SBarry Smith /*------------------------------------------------------------*/
1992d3f70b5SBarry Smith 
2002d3f70b5SBarry Smith /*
2012d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
2022d3f70b5SBarry Smith 
2032d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
2042d3f70b5SBarry Smith */
2054a2ae208SSatish Balay #undef __FUNCT__
2064a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction"
207dfbe8321SBarry Smith PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2082d3f70b5SBarry Smith {
2092d3f70b5SBarry Smith   TS     ts = (TS) ctx;
210ea709b57SSatish Balay   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
211dfbe8321SBarry Smith   PetscErrorCode ierr;
212a7cc72afSBarry Smith   PetscInt i,n;
2132d3f70b5SBarry Smith 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
2152d3f70b5SBarry Smith   /* apply user provided function */
2162d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
2177bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2182d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
2192d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
2202d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
2212d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
2222d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2232d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2242d3f70b5SBarry Smith   }
225888f2ed8SSatish Balay   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
226888f2ed8SSatish Balay   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
227888f2ed8SSatish Balay   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
2282d3f70b5SBarry Smith 
2293a40ed3dSBarry Smith   PetscFunctionReturn(0);
2302d3f70b5SBarry Smith }
2312d3f70b5SBarry Smith 
2322d3f70b5SBarry Smith /*
2332d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2342d3f70b5SBarry Smith 
2352d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2362d3f70b5SBarry Smith */
2374a2ae208SSatish Balay #undef __FUNCT__
2384a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian"
239dfbe8321SBarry Smith PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2402d3f70b5SBarry Smith {
2412d3f70b5SBarry Smith   TS             ts = (TS) ctx;
242dfbe8321SBarry Smith   PetscErrorCode ierr;
2432d3f70b5SBarry Smith 
2443a40ed3dSBarry Smith   PetscFunctionBegin;
2452d3f70b5SBarry Smith   /* construct users Jacobian */
246a7a1495cSBarry Smith   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
2472d3f70b5SBarry Smith 
248614f654bSBarry Smith   /* shift and scale Jacobian */
249*b0f6e734SBarry Smith   ierr = TSScaleShiftMatrices(ts,*AA,*BB,*str);CHKERRQ(ierr);
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
2512d3f70b5SBarry Smith }
2522d3f70b5SBarry Smith 
2532d3f70b5SBarry Smith 
2544a2ae208SSatish Balay #undef __FUNCT__
2554a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2566849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2572d3f70b5SBarry Smith {
2587bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
259dfbe8321SBarry Smith   PetscErrorCode ierr;
2602d3f70b5SBarry Smith 
2613a40ed3dSBarry Smith   PetscFunctionBegin;
262e6e5fe25SBarry Smith   /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */
2637bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2647bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2657bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2667bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
2673a40ed3dSBarry Smith   PetscFunctionReturn(0);
2682d3f70b5SBarry Smith }
2692d3f70b5SBarry Smith /*------------------------------------------------------------*/
2702d3f70b5SBarry Smith 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultMonitor"
273a7cc72afSBarry Smith PetscErrorCode TSPseudoDefaultMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
2742d3f70b5SBarry Smith {
2757bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
276dfbe8321SBarry Smith   PetscErrorCode ierr;
2772d3f70b5SBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionBegin;
27977431f27SBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm,"TS %D dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2812d3f70b5SBarry Smith }
2822d3f70b5SBarry Smith 
2834a2ae208SSatish Balay #undef __FUNCT__
2844a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
2856849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
2862d3f70b5SBarry Smith {
2874bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
288dfbe8321SBarry Smith   PetscErrorCode ierr;
289f1af5d2fSBarry Smith   PetscTruth flg;
2902d3f70b5SBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
2922d3f70b5SBarry Smith 
293b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
294b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
2952d3f70b5SBarry Smith     if (flg) {
296329f5518SBarry Smith       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
29728aa8177SBarry Smith     }
298b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
299ca4b7087SBarry Smith     if (flg) {
300ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
301ca4b7087SBarry Smith     }
302419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
303b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3042d3f70b5SBarry Smith 
3053a40ed3dSBarry Smith   PetscFunctionReturn(0);
3062d3f70b5SBarry Smith }
3072d3f70b5SBarry Smith 
3084a2ae208SSatish Balay #undef __FUNCT__
3094a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
3106849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3112d3f70b5SBarry Smith {
3123a40ed3dSBarry Smith   PetscFunctionBegin;
3133a40ed3dSBarry Smith   PetscFunctionReturn(0);
3142d3f70b5SBarry Smith }
3152d3f70b5SBarry Smith 
31682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3174a2ae208SSatish Balay #undef __FUNCT__
3184a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
319ac226902SBarry Smith /*@C
32082bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
32182bf6240SBarry Smith    last timestep.
32282bf6240SBarry Smith 
32315091d37SBarry Smith    Collective on TS
32415091d37SBarry Smith 
32582bf6240SBarry Smith    Input Parameters:
32615091d37SBarry Smith +  ts - timestep context
32782bf6240SBarry Smith .  dt - user-defined function to verify timestep
32815091d37SBarry Smith -  ctx - [optional] user-defined context for private data
32982bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
33082bf6240SBarry Smith 
33115091d37SBarry Smith    Level: advanced
332fee21e36SBarry Smith 
33382bf6240SBarry Smith    Calling sequence of func:
334a7cc72afSBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscTruth *flag);
33582bf6240SBarry Smith 
33682bf6240SBarry Smith .  update - latest solution vector
33782bf6240SBarry Smith .  ctx - [optional] timestep context
33882bf6240SBarry Smith .  newdt - the timestep to use for the next step
33982bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
34082bf6240SBarry Smith 
34182bf6240SBarry Smith    Notes:
34282bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
34382bf6240SBarry Smith    during the timestepping process.
34482bf6240SBarry Smith 
34582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
34682bf6240SBarry Smith 
34782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
34882bf6240SBarry Smith @*/
349a7cc72afSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscTruth*),void* ctx)
35082bf6240SBarry Smith {
351a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscTruth *),void *);
35282bf6240SBarry Smith 
35382bf6240SBarry Smith   PetscFunctionBegin;
3544482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
35582bf6240SBarry Smith 
356c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
35782bf6240SBarry Smith   if (f) {
35882bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
35982bf6240SBarry Smith   }
36082bf6240SBarry Smith   PetscFunctionReturn(0);
36182bf6240SBarry Smith }
36282bf6240SBarry Smith 
3634a2ae208SSatish Balay #undef __FUNCT__
3644a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
36582bf6240SBarry Smith /*@
36682bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
36782bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
36882bf6240SBarry Smith 
369fee21e36SBarry Smith    Collective on TS
370fee21e36SBarry Smith 
37115091d37SBarry Smith     Input Parameters:
37215091d37SBarry Smith +   ts - the timestep context
37315091d37SBarry Smith -   inc - the scaling factor >= 1.0
37415091d37SBarry Smith 
37582bf6240SBarry Smith     Options Database Key:
37682bf6240SBarry Smith $    -ts_pseudo_increment <increment>
37782bf6240SBarry Smith 
37815091d37SBarry Smith     Level: advanced
37915091d37SBarry Smith 
38082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
38182bf6240SBarry Smith 
38282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
38382bf6240SBarry Smith @*/
384dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
38582bf6240SBarry Smith {
386dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscReal);
38782bf6240SBarry Smith 
38882bf6240SBarry Smith   PetscFunctionBegin;
3894482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
39082bf6240SBarry Smith 
391c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
39282bf6240SBarry Smith   if (f) {
39382bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
39482bf6240SBarry Smith   }
39582bf6240SBarry Smith   PetscFunctionReturn(0);
39682bf6240SBarry Smith }
39782bf6240SBarry Smith 
3984a2ae208SSatish Balay #undef __FUNCT__
3994a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
40082bf6240SBarry Smith /*@
40182bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
40282bf6240SBarry Smith     is computed via the formula
40382bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
40482bf6240SBarry Smith       rather than the default update,
40582bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
40682bf6240SBarry Smith 
40715091d37SBarry Smith    Collective on TS
40815091d37SBarry Smith 
40982bf6240SBarry Smith     Input Parameter:
41082bf6240SBarry Smith .   ts - the timestep context
41182bf6240SBarry Smith 
41282bf6240SBarry Smith     Options Database Key:
41382bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
41482bf6240SBarry Smith 
41515091d37SBarry Smith     Level: advanced
41615091d37SBarry Smith 
41782bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
41882bf6240SBarry Smith 
41982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42082bf6240SBarry Smith @*/
421dfbe8321SBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
42282bf6240SBarry Smith {
423dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS);
42482bf6240SBarry Smith 
42582bf6240SBarry Smith   PetscFunctionBegin;
4264482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
42782bf6240SBarry Smith 
428c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
42982bf6240SBarry Smith   if (f) {
43082bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
43182bf6240SBarry Smith   }
43282bf6240SBarry Smith   PetscFunctionReturn(0);
43382bf6240SBarry Smith }
43482bf6240SBarry Smith 
43582bf6240SBarry Smith 
4364a2ae208SSatish Balay #undef __FUNCT__
4374a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
438ac226902SBarry Smith /*@C
43982bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
44082bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
44182bf6240SBarry Smith 
44215091d37SBarry Smith    Collective on TS
44315091d37SBarry Smith 
44482bf6240SBarry Smith    Input Parameters:
44515091d37SBarry Smith +  ts - timestep context
44682bf6240SBarry Smith .  dt - function to compute timestep
44715091d37SBarry Smith -  ctx - [optional] user-defined context for private data
44882bf6240SBarry Smith          required by the function (may be PETSC_NULL)
44982bf6240SBarry Smith 
45015091d37SBarry Smith    Level: intermediate
451fee21e36SBarry Smith 
45282bf6240SBarry Smith    Calling sequence of func:
45387828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
45482bf6240SBarry Smith 
45582bf6240SBarry Smith .  newdt - the newly computed timestep
45682bf6240SBarry Smith .  ctx - [optional] timestep context
45782bf6240SBarry Smith 
45882bf6240SBarry Smith    Notes:
45982bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
46082bf6240SBarry Smith    during the timestepping process.
46182bf6240SBarry Smith 
46282bf6240SBarry Smith .keywords: timestep, pseudo, set
46382bf6240SBarry Smith 
46482bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
46582bf6240SBarry Smith @*/
4666849ba73SBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
46782bf6240SBarry Smith {
4686849ba73SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *);
46982bf6240SBarry Smith 
47082bf6240SBarry Smith   PetscFunctionBegin;
4714482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
47282bf6240SBarry Smith 
473c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
47482bf6240SBarry Smith   if (f) {
47582bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
47682bf6240SBarry Smith   }
47782bf6240SBarry Smith   PetscFunctionReturn(0);
47882bf6240SBarry Smith }
47982bf6240SBarry Smith 
48082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
48182bf6240SBarry Smith 
482a7cc72afSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
483fb2e594dSBarry Smith EXTERN_C_BEGIN
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
486dfbe8321SBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
48782bf6240SBarry Smith {
48882bf6240SBarry Smith   TS_Pseudo *pseudo;
48982bf6240SBarry Smith 
49082bf6240SBarry Smith   PetscFunctionBegin;
49182bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
49282bf6240SBarry Smith   pseudo->verify      = dt;
49382bf6240SBarry Smith   pseudo->verifyctx   = ctx;
49482bf6240SBarry Smith   PetscFunctionReturn(0);
49582bf6240SBarry Smith }
496fb2e594dSBarry Smith EXTERN_C_END
49782bf6240SBarry Smith 
498fb2e594dSBarry Smith EXTERN_C_BEGIN
4994a2ae208SSatish Balay #undef __FUNCT__
5004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
501dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
50282bf6240SBarry Smith {
5034bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
50482bf6240SBarry Smith 
50582bf6240SBarry Smith   PetscFunctionBegin;
50682bf6240SBarry Smith   pseudo->dt_increment = inc;
50782bf6240SBarry Smith   PetscFunctionReturn(0);
50882bf6240SBarry Smith }
509fb2e594dSBarry Smith EXTERN_C_END
51082bf6240SBarry Smith 
511fb2e594dSBarry Smith EXTERN_C_BEGIN
5124a2ae208SSatish Balay #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
514dfbe8321SBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
51582bf6240SBarry Smith {
5164bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
51782bf6240SBarry Smith 
51882bf6240SBarry Smith   PetscFunctionBegin;
5194bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
52082bf6240SBarry Smith   PetscFunctionReturn(0);
52182bf6240SBarry Smith }
522fb2e594dSBarry Smith EXTERN_C_END
52382bf6240SBarry Smith 
5246849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
525fb2e594dSBarry Smith EXTERN_C_BEGIN
5264a2ae208SSatish Balay #undef __FUNCT__
5274a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
528dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
52982bf6240SBarry Smith {
5304bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53182bf6240SBarry Smith 
53282bf6240SBarry Smith   PetscFunctionBegin;
53382bf6240SBarry Smith   pseudo->dt      = dt;
53482bf6240SBarry Smith   pseudo->dtctx   = ctx;
53582bf6240SBarry Smith   PetscFunctionReturn(0);
53682bf6240SBarry Smith }
537fb2e594dSBarry Smith EXTERN_C_END
53882bf6240SBarry Smith 
53982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
54082bf6240SBarry Smith 
541fb2e594dSBarry Smith EXTERN_C_BEGIN
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
544dfbe8321SBarry Smith PetscErrorCode TSCreate_Pseudo(TS ts)
5452d3f70b5SBarry Smith {
5467bf11e45SBarry Smith   TS_Pseudo  *pseudo;
547dfbe8321SBarry Smith   PetscErrorCode ierr;
5482d3f70b5SBarry Smith 
5493a40ed3dSBarry Smith   PetscFunctionBegin;
550000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
551000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
5522d3f70b5SBarry Smith 
5532d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
55429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
5552d3f70b5SBarry Smith   }
5562d3f70b5SBarry Smith   if (!ts->A) {
55729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
5582d3f70b5SBarry Smith   }
559273d9f13SBarry Smith 
560000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
561000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
562000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
5637bf11e45SBarry Smith 
5647bf11e45SBarry Smith   /* create the required nonlinear solver context */
5654b27c08aSLois Curfman McInnes   ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr);
5662d3f70b5SBarry Smith 
567b0a32e0cSBarry Smith   ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr);
568b0a32e0cSBarry Smith   PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
569eed86810SBarry Smith 
570549d3d68SSatish Balay   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
5717bf11e45SBarry Smith   ts->data = (void*)pseudo;
5722d3f70b5SBarry Smith 
57328aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
5744bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
57528aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
57682bf6240SBarry Smith 
577f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
578e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
5790c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
580f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
581e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
5820c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
583f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
584e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
5850c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
5860c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
5870c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
5880c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
5902d3f70b5SBarry Smith }
591fb2e594dSBarry Smith EXTERN_C_END
5922d3f70b5SBarry Smith 
5934a2ae208SSatish Balay #undef __FUNCT__
5944a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
59582bf6240SBarry Smith /*@C
59682bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
59782bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
59828aa8177SBarry Smith 
59915091d37SBarry Smith    Collective on TS
60015091d37SBarry Smith 
60128aa8177SBarry Smith    Input Parameters:
60228aa8177SBarry Smith .  ts - the timestep context
60382bf6240SBarry Smith .  dtctx - unused timestep context
60428aa8177SBarry Smith 
60582bf6240SBarry Smith    Output Parameter:
60682bf6240SBarry Smith .  newdt - the timestep to use for the next step
60728aa8177SBarry Smith 
60815091d37SBarry Smith    Level: advanced
60915091d37SBarry Smith 
61082bf6240SBarry Smith .keywords: timestep, pseudo, default
611564e8f4eSLois Curfman McInnes 
61282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
61328aa8177SBarry Smith @*/
614dfbe8321SBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
61528aa8177SBarry Smith {
61682bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61787828ca2SBarry Smith   PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
618dfbe8321SBarry Smith   PetscErrorCode ierr;
61928aa8177SBarry Smith 
6203a40ed3dSBarry Smith   PetscFunctionBegin;
62182bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
62282bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
62382bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
62482bf6240SBarry Smith     /* first time through so compute initial function norm */
62582bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
62682bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
62782bf6240SBarry Smith   }
62882bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
62982bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
630004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
63182bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
63282bf6240SBarry Smith   } else {
63382bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
63482bf6240SBarry Smith   }
63582bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6363a40ed3dSBarry Smith   PetscFunctionReturn(0);
63728aa8177SBarry Smith }
6382d3f70b5SBarry Smith 
639004884caSBarry Smith 
640004884caSBarry Smith 
641004884caSBarry Smith 
642004884caSBarry Smith 
643004884caSBarry Smith 
644004884caSBarry Smith 
645004884caSBarry Smith 
646004884caSBarry Smith 
647004884caSBarry Smith 
648004884caSBarry Smith 
649004884caSBarry Smith 
650004884caSBarry Smith 
651004884caSBarry Smith 
652004884caSBarry Smith 
653004884caSBarry Smith 
654004884caSBarry Smith 
655