xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
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 
13*6849ba73SBarry Smith   PetscErrorCode (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
142d3f70b5SBarry Smith   void   *dtctx;
15*6849ba73SBarry Smith   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,int*); /* 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 @*/
91dfbe8321SBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag)
927bf11e45SBarry Smith {
933a40ed3dSBarry Smith   PetscFunctionBegin;
947bf11e45SBarry Smith   *flag = 1;
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 @*/
124dfbe8321SBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag)
1257bf11e45SBarry Smith {
1267bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
127dfbe8321SBarry Smith   PetscErrorCode ierr;
1287bf11e45SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionBegin;
1303a40ed3dSBarry Smith   if (!pseudo->verify) {*flag = 1; 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"
141*6849ba73SBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime)
1422d3f70b5SBarry Smith {
1432d3f70b5SBarry Smith   Vec       sol = ts->vec_sol;
144dfbe8321SBarry Smith   PetscErrorCode ierr;
145dfbe8321SBarry Smith   int i,max_steps = ts->max_steps,its,ok,lits;
1467bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
14787828ca2SBarry Smith   PetscReal current_time_step;
1482d3f70b5SBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
1502d3f70b5SBarry Smith   *steps = -ts->steps;
1512d3f70b5SBarry Smith 
1527bf11e45SBarry Smith   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
1532d3f70b5SBarry Smith   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
1547bf11e45SBarry Smith     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
155e6e5fe25SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1567bf11e45SBarry Smith     current_time_step = ts->time_step;
157e6e5fe25SBarry Smith     while (PETSC_TRUE) {
1587bf11e45SBarry Smith       ts->ptime  += current_time_step;
159c9780b6fSBarry Smith       ierr = SNESSolve(ts->snes,pseudo->update);CHKERRQ(ierr);
160f2267985SLois Curfman McInnes       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
1619f954729SBarry Smith       ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
162c9780b6fSBarry Smith       ts->nonlinear_its += its; ts->linear_its += lits;
1637bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
1647bf11e45SBarry Smith       if (ok) break;
1657bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1667bf11e45SBarry Smith       current_time_step = ts->time_step;
1677bf11e45SBarry Smith     }
1687bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
1692d3f70b5SBarry Smith     ts->steps++;
1702d3f70b5SBarry Smith   }
171e6e5fe25SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
172e6e5fe25SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
173e6e5fe25SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1742d3f70b5SBarry Smith 
1752d3f70b5SBarry Smith   *steps += ts->steps;
176142b95e3SSatish Balay   *ptime  = ts->ptime;
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1782d3f70b5SBarry Smith }
1792d3f70b5SBarry Smith 
1802d3f70b5SBarry Smith /*------------------------------------------------------------*/
1814a2ae208SSatish Balay #undef __FUNCT__
1824a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo"
183*6849ba73SBarry Smith static PetscErrorCode TSDestroy_Pseudo(TS ts)
1842d3f70b5SBarry Smith {
1857bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
186dfbe8321SBarry Smith   PetscErrorCode ierr;
1872d3f70b5SBarry Smith 
1883a40ed3dSBarry Smith   PetscFunctionBegin;
189e4ef1970SSatish Balay   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
1907bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
1917bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
192606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1933a40ed3dSBarry Smith   PetscFunctionReturn(0);
1942d3f70b5SBarry Smith }
1952d3f70b5SBarry Smith 
1962d3f70b5SBarry Smith 
1972d3f70b5SBarry Smith /*------------------------------------------------------------*/
1982d3f70b5SBarry Smith 
1992d3f70b5SBarry Smith /*
2002d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
2012d3f70b5SBarry Smith 
2022d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
2032d3f70b5SBarry Smith */
2044a2ae208SSatish Balay #undef __FUNCT__
2054a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction"
206dfbe8321SBarry Smith PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2072d3f70b5SBarry Smith {
2082d3f70b5SBarry Smith   TS     ts = (TS) ctx;
209ea709b57SSatish Balay   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
210dfbe8321SBarry Smith   PetscErrorCode ierr;
211dfbe8321SBarry Smith   int i,n;
2122d3f70b5SBarry Smith 
2133a40ed3dSBarry Smith   PetscFunctionBegin;
2142d3f70b5SBarry Smith   /* apply user provided function */
2152d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
2167bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2172d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
2182d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
2192d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
2202d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
2212d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2222d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2232d3f70b5SBarry Smith   }
224888f2ed8SSatish Balay   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
225888f2ed8SSatish Balay   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
226888f2ed8SSatish Balay   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
2272d3f70b5SBarry Smith 
2283a40ed3dSBarry Smith   PetscFunctionReturn(0);
2292d3f70b5SBarry Smith }
2302d3f70b5SBarry Smith 
2312d3f70b5SBarry Smith /*
2322d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2332d3f70b5SBarry Smith 
2342d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2352d3f70b5SBarry Smith */
2364a2ae208SSatish Balay #undef __FUNCT__
2374a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian"
238dfbe8321SBarry Smith PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2392d3f70b5SBarry Smith {
2402d3f70b5SBarry Smith   TS            ts = (TS) ctx;
241dfbe8321SBarry Smith   PetscErrorCode ierr;
242ea709b57SSatish Balay   PetscScalar   mone = -1.0,mdt = 1.0/ts->time_step;
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 */
2492d3f70b5SBarry Smith   ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
2502d3f70b5SBarry Smith   ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
251614f654bSBarry Smith   if (*BB != *AA && *str != SAME_PRECONDITIONER) {
2522d3f70b5SBarry Smith     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
2532d3f70b5SBarry Smith     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
2542d3f70b5SBarry Smith   }
2552d3f70b5SBarry Smith 
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2572d3f70b5SBarry Smith }
2582d3f70b5SBarry Smith 
2592d3f70b5SBarry Smith 
2604a2ae208SSatish Balay #undef __FUNCT__
2614a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
262*6849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts)
2632d3f70b5SBarry Smith {
2647bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
265dfbe8321SBarry Smith   PetscErrorCode ierr;
2662d3f70b5SBarry Smith 
2673a40ed3dSBarry Smith   PetscFunctionBegin;
268e6e5fe25SBarry Smith   /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */
2697bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2707bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2717bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2727bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2742d3f70b5SBarry Smith }
2752d3f70b5SBarry Smith /*------------------------------------------------------------*/
2762d3f70b5SBarry Smith 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultMonitor"
279dfbe8321SBarry Smith PetscErrorCode TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
2802d3f70b5SBarry Smith {
2817bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
2832d3f70b5SBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
285142b95e3SSatish Balay   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
2863a40ed3dSBarry Smith   PetscFunctionReturn(0);
2872d3f70b5SBarry Smith }
2882d3f70b5SBarry Smith 
2894a2ae208SSatish Balay #undef __FUNCT__
2904a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
291*6849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
2922d3f70b5SBarry Smith {
2934bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
294dfbe8321SBarry Smith   PetscErrorCode ierr;
295f1af5d2fSBarry Smith   PetscTruth flg;
2962d3f70b5SBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2982d3f70b5SBarry Smith 
299b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
300b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
3012d3f70b5SBarry Smith     if (flg) {
302329f5518SBarry Smith       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
30328aa8177SBarry Smith     }
304b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
305ca4b7087SBarry Smith     if (flg) {
306ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
307ca4b7087SBarry Smith     }
308419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
309b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3102d3f70b5SBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
3122d3f70b5SBarry Smith }
3132d3f70b5SBarry Smith 
3144a2ae208SSatish Balay #undef __FUNCT__
3154a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
316*6849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
3172d3f70b5SBarry Smith {
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
3202d3f70b5SBarry Smith }
3212d3f70b5SBarry Smith 
32282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3234a2ae208SSatish Balay #undef __FUNCT__
3244a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
325ac226902SBarry Smith /*@C
32682bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
32782bf6240SBarry Smith    last timestep.
32882bf6240SBarry Smith 
32915091d37SBarry Smith    Collective on TS
33015091d37SBarry Smith 
33182bf6240SBarry Smith    Input Parameters:
33215091d37SBarry Smith +  ts - timestep context
33382bf6240SBarry Smith .  dt - user-defined function to verify timestep
33415091d37SBarry Smith -  ctx - [optional] user-defined context for private data
33582bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
33682bf6240SBarry Smith 
33715091d37SBarry Smith    Level: advanced
338fee21e36SBarry Smith 
33982bf6240SBarry Smith    Calling sequence of func:
34087828ca2SBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
34182bf6240SBarry Smith 
34282bf6240SBarry Smith .  update - latest solution vector
34382bf6240SBarry Smith .  ctx - [optional] timestep context
34482bf6240SBarry Smith .  newdt - the timestep to use for the next step
34582bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
34682bf6240SBarry Smith 
34782bf6240SBarry Smith    Notes:
34882bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
34982bf6240SBarry Smith    during the timestepping process.
35082bf6240SBarry Smith 
35182bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
35282bf6240SBarry Smith 
35382bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
35482bf6240SBarry Smith @*/
355*6849ba73SBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
35682bf6240SBarry Smith {
357*6849ba73SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,int *),void *);
35882bf6240SBarry Smith 
35982bf6240SBarry Smith   PetscFunctionBegin;
3604482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
36182bf6240SBarry Smith 
362c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
36382bf6240SBarry Smith   if (f) {
36482bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
36582bf6240SBarry Smith   }
36682bf6240SBarry Smith   PetscFunctionReturn(0);
36782bf6240SBarry Smith }
36882bf6240SBarry Smith 
3694a2ae208SSatish Balay #undef __FUNCT__
3704a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
37182bf6240SBarry Smith /*@
37282bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
37382bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
37482bf6240SBarry Smith 
375fee21e36SBarry Smith    Collective on TS
376fee21e36SBarry Smith 
37715091d37SBarry Smith     Input Parameters:
37815091d37SBarry Smith +   ts - the timestep context
37915091d37SBarry Smith -   inc - the scaling factor >= 1.0
38015091d37SBarry Smith 
38182bf6240SBarry Smith     Options Database Key:
38282bf6240SBarry Smith $    -ts_pseudo_increment <increment>
38382bf6240SBarry Smith 
38415091d37SBarry Smith     Level: advanced
38515091d37SBarry Smith 
38682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
38782bf6240SBarry Smith 
38882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
38982bf6240SBarry Smith @*/
390dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
39182bf6240SBarry Smith {
392dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscReal);
39382bf6240SBarry Smith 
39482bf6240SBarry Smith   PetscFunctionBegin;
3954482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
39682bf6240SBarry Smith 
397c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
39882bf6240SBarry Smith   if (f) {
39982bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
40082bf6240SBarry Smith   }
40182bf6240SBarry Smith   PetscFunctionReturn(0);
40282bf6240SBarry Smith }
40382bf6240SBarry Smith 
4044a2ae208SSatish Balay #undef __FUNCT__
4054a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
40682bf6240SBarry Smith /*@
40782bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
40882bf6240SBarry Smith     is computed via the formula
40982bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
41082bf6240SBarry Smith       rather than the default update,
41182bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
41282bf6240SBarry Smith 
41315091d37SBarry Smith    Collective on TS
41415091d37SBarry Smith 
41582bf6240SBarry Smith     Input Parameter:
41682bf6240SBarry Smith .   ts - the timestep context
41782bf6240SBarry Smith 
41882bf6240SBarry Smith     Options Database Key:
41982bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
42082bf6240SBarry Smith 
42115091d37SBarry Smith     Level: advanced
42215091d37SBarry Smith 
42382bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
42482bf6240SBarry Smith 
42582bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42682bf6240SBarry Smith @*/
427dfbe8321SBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
42882bf6240SBarry Smith {
429dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(TS);
43082bf6240SBarry Smith 
43182bf6240SBarry Smith   PetscFunctionBegin;
4324482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
43382bf6240SBarry Smith 
434c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
43582bf6240SBarry Smith   if (f) {
43682bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
43782bf6240SBarry Smith   }
43882bf6240SBarry Smith   PetscFunctionReturn(0);
43982bf6240SBarry Smith }
44082bf6240SBarry Smith 
44182bf6240SBarry Smith 
4424a2ae208SSatish Balay #undef __FUNCT__
4434a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
444ac226902SBarry Smith /*@C
44582bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
44682bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
44782bf6240SBarry Smith 
44815091d37SBarry Smith    Collective on TS
44915091d37SBarry Smith 
45082bf6240SBarry Smith    Input Parameters:
45115091d37SBarry Smith +  ts - timestep context
45282bf6240SBarry Smith .  dt - function to compute timestep
45315091d37SBarry Smith -  ctx - [optional] user-defined context for private data
45482bf6240SBarry Smith          required by the function (may be PETSC_NULL)
45582bf6240SBarry Smith 
45615091d37SBarry Smith    Level: intermediate
457fee21e36SBarry Smith 
45882bf6240SBarry Smith    Calling sequence of func:
45987828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
46082bf6240SBarry Smith 
46182bf6240SBarry Smith .  newdt - the newly computed timestep
46282bf6240SBarry Smith .  ctx - [optional] timestep context
46382bf6240SBarry Smith 
46482bf6240SBarry Smith    Notes:
46582bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
46682bf6240SBarry Smith    during the timestepping process.
46782bf6240SBarry Smith 
46882bf6240SBarry Smith .keywords: timestep, pseudo, set
46982bf6240SBarry Smith 
47082bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
47182bf6240SBarry Smith @*/
472*6849ba73SBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
47382bf6240SBarry Smith {
474*6849ba73SBarry Smith   PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *);
47582bf6240SBarry Smith 
47682bf6240SBarry Smith   PetscFunctionBegin;
4774482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
47882bf6240SBarry Smith 
479c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
48082bf6240SBarry Smith   if (f) {
48182bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
48282bf6240SBarry Smith   }
48382bf6240SBarry Smith   PetscFunctionReturn(0);
48482bf6240SBarry Smith }
48582bf6240SBarry Smith 
48682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
48782bf6240SBarry Smith 
488*6849ba73SBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,int*); /* force argument to next function to not be extern C*/
489fb2e594dSBarry Smith EXTERN_C_BEGIN
4904a2ae208SSatish Balay #undef __FUNCT__
4914a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
492dfbe8321SBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
49382bf6240SBarry Smith {
49482bf6240SBarry Smith   TS_Pseudo *pseudo;
49582bf6240SBarry Smith 
49682bf6240SBarry Smith   PetscFunctionBegin;
49782bf6240SBarry Smith   pseudo              = (TS_Pseudo*)ts->data;
49882bf6240SBarry Smith   pseudo->verify      = dt;
49982bf6240SBarry Smith   pseudo->verifyctx   = ctx;
50082bf6240SBarry Smith   PetscFunctionReturn(0);
50182bf6240SBarry Smith }
502fb2e594dSBarry Smith EXTERN_C_END
50382bf6240SBarry Smith 
504fb2e594dSBarry Smith EXTERN_C_BEGIN
5054a2ae208SSatish Balay #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
507dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
50882bf6240SBarry Smith {
5094bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
51082bf6240SBarry Smith 
51182bf6240SBarry Smith   PetscFunctionBegin;
51282bf6240SBarry Smith   pseudo->dt_increment = inc;
51382bf6240SBarry Smith   PetscFunctionReturn(0);
51482bf6240SBarry Smith }
515fb2e594dSBarry Smith EXTERN_C_END
51682bf6240SBarry Smith 
517fb2e594dSBarry Smith EXTERN_C_BEGIN
5184a2ae208SSatish Balay #undef __FUNCT__
5194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
520dfbe8321SBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
52182bf6240SBarry Smith {
5224bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
52382bf6240SBarry Smith 
52482bf6240SBarry Smith   PetscFunctionBegin;
5254bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
52682bf6240SBarry Smith   PetscFunctionReturn(0);
52782bf6240SBarry Smith }
528fb2e594dSBarry Smith EXTERN_C_END
52982bf6240SBarry Smith 
530*6849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
531fb2e594dSBarry Smith EXTERN_C_BEGIN
5324a2ae208SSatish Balay #undef __FUNCT__
5334a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
534dfbe8321SBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
53582bf6240SBarry Smith {
5364bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53782bf6240SBarry Smith 
53882bf6240SBarry Smith   PetscFunctionBegin;
53982bf6240SBarry Smith   pseudo->dt      = dt;
54082bf6240SBarry Smith   pseudo->dtctx   = ctx;
54182bf6240SBarry Smith   PetscFunctionReturn(0);
54282bf6240SBarry Smith }
543fb2e594dSBarry Smith EXTERN_C_END
54482bf6240SBarry Smith 
54582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
54682bf6240SBarry Smith 
547fb2e594dSBarry Smith EXTERN_C_BEGIN
5484a2ae208SSatish Balay #undef __FUNCT__
5494a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
550dfbe8321SBarry Smith PetscErrorCode TSCreate_Pseudo(TS ts)
5512d3f70b5SBarry Smith {
5527bf11e45SBarry Smith   TS_Pseudo  *pseudo;
553dfbe8321SBarry Smith   PetscErrorCode ierr;
5542d3f70b5SBarry Smith 
5553a40ed3dSBarry Smith   PetscFunctionBegin;
556000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Pseudo;
557000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Pseudo;
5582d3f70b5SBarry Smith 
5592d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
56029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
5612d3f70b5SBarry Smith   }
5622d3f70b5SBarry Smith   if (!ts->A) {
56329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
5642d3f70b5SBarry Smith   }
565273d9f13SBarry Smith 
566000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Pseudo;
567000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Pseudo;
568000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
5697bf11e45SBarry Smith 
5707bf11e45SBarry Smith   /* create the required nonlinear solver context */
5714b27c08aSLois Curfman McInnes   ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr);
5722d3f70b5SBarry Smith 
573b0a32e0cSBarry Smith   ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr);
574b0a32e0cSBarry Smith   PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
575eed86810SBarry Smith 
576549d3d68SSatish Balay   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
5777bf11e45SBarry Smith   ts->data = (void*)pseudo;
5782d3f70b5SBarry Smith 
57928aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
5804bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
58128aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
58282bf6240SBarry Smith 
583f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
584e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
5850c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
586f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
587e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
5880c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
589f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
590e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
5910c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
5920c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
5930c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
5940c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
5953a40ed3dSBarry Smith   PetscFunctionReturn(0);
5962d3f70b5SBarry Smith }
597fb2e594dSBarry Smith EXTERN_C_END
5982d3f70b5SBarry Smith 
5994a2ae208SSatish Balay #undef __FUNCT__
6004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
60182bf6240SBarry Smith /*@C
60282bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
60382bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
60428aa8177SBarry Smith 
60515091d37SBarry Smith    Collective on TS
60615091d37SBarry Smith 
60728aa8177SBarry Smith    Input Parameters:
60828aa8177SBarry Smith .  ts - the timestep context
60982bf6240SBarry Smith .  dtctx - unused timestep context
61028aa8177SBarry Smith 
61182bf6240SBarry Smith    Output Parameter:
61282bf6240SBarry Smith .  newdt - the timestep to use for the next step
61328aa8177SBarry Smith 
61415091d37SBarry Smith    Level: advanced
61515091d37SBarry Smith 
61682bf6240SBarry Smith .keywords: timestep, pseudo, default
617564e8f4eSLois Curfman McInnes 
61882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
61928aa8177SBarry Smith @*/
620dfbe8321SBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
62128aa8177SBarry Smith {
62282bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
62387828ca2SBarry Smith   PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
624dfbe8321SBarry Smith   PetscErrorCode ierr;
62528aa8177SBarry Smith 
6263a40ed3dSBarry Smith   PetscFunctionBegin;
62782bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
62882bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
62982bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
63082bf6240SBarry Smith     /* first time through so compute initial function norm */
63182bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
63282bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
63382bf6240SBarry Smith   }
63482bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
63582bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
636004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
63782bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
63882bf6240SBarry Smith   } else {
63982bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
64082bf6240SBarry Smith   }
64182bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6423a40ed3dSBarry Smith   PetscFunctionReturn(0);
64328aa8177SBarry Smith }
6442d3f70b5SBarry Smith 
645004884caSBarry Smith 
646004884caSBarry Smith 
647004884caSBarry Smith 
648004884caSBarry Smith 
649004884caSBarry Smith 
650004884caSBarry Smith 
651004884caSBarry Smith 
652004884caSBarry Smith 
653004884caSBarry Smith 
654004884caSBarry Smith 
655004884caSBarry Smith 
656004884caSBarry Smith 
657004884caSBarry Smith 
658004884caSBarry Smith 
659004884caSBarry Smith 
660004884caSBarry Smith 
661