xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 614f654b7d87221fea7c4c05971e5a5b2bfa5406)
1*614f654bSBarry Smith /*$Id: posindep.c,v 1.58 2001/08/07 17:20:38 balay Exp bsmith $*/
22d3f70b5SBarry Smith /*
3fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
42d3f70b5SBarry Smith */
5e090d566SSatish Balay #include "src/ts/tsimpl.h"                /*I   "petscts.h"   I*/
62d3f70b5SBarry Smith 
72d3f70b5SBarry Smith typedef struct {
82d3f70b5SBarry Smith   Vec  update;      /* work vector where new solution is formed */
92d3f70b5SBarry Smith   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
102d3f70b5SBarry Smith   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
112d3f70b5SBarry Smith 
122d3f70b5SBarry Smith   /* information used for Pseudo-timestepping */
132d3f70b5SBarry Smith 
1487828ca2SBarry Smith   int    (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
152d3f70b5SBarry Smith   void   *dtctx;
1687828ca2SBarry Smith   int    (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */
177bf11e45SBarry Smith   void   *verifyctx;
182d3f70b5SBarry Smith 
1987828ca2SBarry Smith   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
2087828ca2SBarry Smith   PetscReal  fnorm_previous;
2128aa8177SBarry Smith 
2287828ca2SBarry Smith   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
234bbc92c1SBarry Smith   PetscTruth increment_dt_from_initial_dt;
247bf11e45SBarry Smith } TS_Pseudo;
252d3f70b5SBarry Smith 
262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
272d3f70b5SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep"
307bf11e45SBarry Smith /*@
317bf11e45SBarry Smith     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32564e8f4eSLois Curfman McInnes     pseudo-timestepping process.
332d3f70b5SBarry Smith 
3415091d37SBarry Smith     Collective on TS
3515091d37SBarry Smith 
367bf11e45SBarry Smith     Input Parameter:
377bf11e45SBarry Smith .   ts - timestep context
387bf11e45SBarry Smith 
397bf11e45SBarry Smith     Output Parameter:
40fb4a63b6SLois Curfman McInnes .   dt - newly computed timestep
41fb4a63b6SLois Curfman McInnes 
4215091d37SBarry Smith     Level: advanced
43564e8f4eSLois Curfman McInnes 
44564e8f4eSLois Curfman McInnes     Notes:
45564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
46564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetTimeStep().
47564e8f4eSLois Curfman McInnes 
48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute
49564e8f4eSLois Curfman McInnes 
50564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
517bf11e45SBarry Smith @*/
5287828ca2SBarry Smith int TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
537bf11e45SBarry Smith {
547bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
557bf11e45SBarry Smith   int       ierr;
567bf11e45SBarry Smith 
573a40ed3dSBarry Smith   PetscFunctionBegin;
58b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
597bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
60b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
613a40ed3dSBarry Smith   PetscFunctionReturn(0);
627bf11e45SBarry Smith }
637bf11e45SBarry Smith 
647bf11e45SBarry Smith 
657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
687bf11e45SBarry Smith /*@C
69639f9d9dSBarry Smith    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
707bf11e45SBarry Smith 
7115091d37SBarry Smith    Collective on TS
7215091d37SBarry Smith 
737bf11e45SBarry Smith    Input Parameters:
7415091d37SBarry Smith +  ts - the timestep context
757bf11e45SBarry Smith .  dtctx - unused timestep context
7615091d37SBarry Smith -  update - latest solution vector
777bf11e45SBarry Smith 
78564e8f4eSLois Curfman McInnes    Output Parameters:
7915091d37SBarry Smith +  newdt - the timestep to use for the next step
8015091d37SBarry Smith -  flag - flag indicating whether the last time step was acceptable
817bf11e45SBarry Smith 
8215091d37SBarry Smith    Level: advanced
83fee21e36SBarry Smith 
84564e8f4eSLois Curfman McInnes    Note:
85564e8f4eSLois Curfman McInnes    This routine always returns a flag of 1, indicating an acceptable
86564e8f4eSLois Curfman McInnes    timestep.
87564e8f4eSLois Curfman McInnes 
88564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify
89564e8f4eSLois Curfman McInnes 
90564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
917bf11e45SBarry Smith @*/
9287828ca2SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag)
937bf11e45SBarry Smith {
943a40ed3dSBarry Smith   PetscFunctionBegin;
957bf11e45SBarry Smith   *flag = 1;
963a40ed3dSBarry Smith   PetscFunctionReturn(0);
977bf11e45SBarry Smith }
987bf11e45SBarry Smith 
997bf11e45SBarry Smith 
1004a2ae208SSatish Balay #undef __FUNCT__
1014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep"
1027bf11e45SBarry Smith /*@
103564e8f4eSLois Curfman McInnes     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
1047bf11e45SBarry Smith 
10515091d37SBarry Smith     Collective on TS
10615091d37SBarry Smith 
107fb4a63b6SLois Curfman McInnes     Input Parameters:
10815091d37SBarry Smith +   ts - timestep context
10915091d37SBarry Smith -   update - latest solution vector
1107bf11e45SBarry Smith 
111fb4a63b6SLois Curfman McInnes     Output Parameters:
11215091d37SBarry Smith +   dt - newly computed timestep (if it had to shrink)
11315091d37SBarry Smith -   flag - indicates if current timestep was ok
1147bf11e45SBarry Smith 
11515091d37SBarry Smith     Level: advanced
116fee21e36SBarry Smith 
117564e8f4eSLois Curfman McInnes     Notes:
118564e8f4eSLois Curfman McInnes     The routine to be called here to compute the timestep should be
119564e8f4eSLois Curfman McInnes     set by calling TSPseudoSetVerifyTimeStep().
120564e8f4eSLois Curfman McInnes 
121564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify
122564e8f4eSLois Curfman McInnes 
123564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
1247bf11e45SBarry Smith @*/
12587828ca2SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag)
1267bf11e45SBarry Smith {
1277bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
1287bf11e45SBarry Smith   int       ierr;
1297bf11e45SBarry Smith 
1303a40ed3dSBarry Smith   PetscFunctionBegin;
1313a40ed3dSBarry Smith   if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);}
1327bf11e45SBarry Smith 
1337bf11e45SBarry Smith   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
1347bf11e45SBarry Smith 
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1367bf11e45SBarry Smith }
1377bf11e45SBarry Smith 
1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/
1397bf11e45SBarry Smith 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo"
14287828ca2SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime)
1432d3f70b5SBarry Smith {
1442d3f70b5SBarry Smith   Vec       sol = ts->vec_sol;
145f2267985SLois Curfman McInnes   int       ierr,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);
1557bf11e45SBarry Smith     current_time_step = ts->time_step;
1567bf11e45SBarry Smith     while (1) {
1577bf11e45SBarry Smith       ts->ptime  += current_time_step;
1587bf11e45SBarry Smith       ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr);
159f2267985SLois Curfman McInnes       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
160f2267985SLois Curfman McInnes       ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
1617bf11e45SBarry Smith       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
1627bf11e45SBarry Smith       if (ok) break;
1637bf11e45SBarry Smith       ts->ptime        -= current_time_step;
1647bf11e45SBarry Smith       current_time_step = ts->time_step;
1657bf11e45SBarry Smith     }
1667bf11e45SBarry Smith     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
1672d3f70b5SBarry Smith     ts->steps++;
1682d3f70b5SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1692d3f70b5SBarry Smith   }
1702d3f70b5SBarry Smith 
1712d3f70b5SBarry Smith   *steps += ts->steps;
172142b95e3SSatish Balay   *ptime  = ts->ptime;
1733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1742d3f70b5SBarry Smith }
1752d3f70b5SBarry Smith 
1762d3f70b5SBarry Smith /*------------------------------------------------------------*/
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Pseudo"
179e1311b90SBarry Smith static int TSDestroy_Pseudo(TS ts)
1802d3f70b5SBarry Smith {
1817bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
1822d3f70b5SBarry Smith   int       ierr;
1832d3f70b5SBarry Smith 
1843a40ed3dSBarry Smith   PetscFunctionBegin;
185e4ef1970SSatish Balay   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
1867bf11e45SBarry Smith   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
1877bf11e45SBarry Smith   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
188606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1893a40ed3dSBarry Smith   PetscFunctionReturn(0);
1902d3f70b5SBarry Smith }
1912d3f70b5SBarry Smith 
1922d3f70b5SBarry Smith 
1932d3f70b5SBarry Smith /*------------------------------------------------------------*/
1942d3f70b5SBarry Smith 
1952d3f70b5SBarry Smith /*
1962d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
1972d3f70b5SBarry Smith 
1982d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
1992d3f70b5SBarry Smith */
2004a2ae208SSatish Balay #undef __FUNCT__
2014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoFunction"
2027bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2032d3f70b5SBarry Smith {
2042d3f70b5SBarry Smith   TS     ts = (TS) ctx;
205ea709b57SSatish Balay   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
2062d3f70b5SBarry Smith   int    ierr,i,n;
2072d3f70b5SBarry Smith 
2083a40ed3dSBarry Smith   PetscFunctionBegin;
2092d3f70b5SBarry Smith   /* apply user provided function */
2102d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
2117bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2122d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
2132d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
2142d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
2152d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
2162d3f70b5SBarry Smith   for (i=0; i<n; i++) {
2172d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2182d3f70b5SBarry Smith   }
219888f2ed8SSatish Balay   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
220888f2ed8SSatish Balay   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
221888f2ed8SSatish Balay   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
2222d3f70b5SBarry Smith 
2233a40ed3dSBarry Smith   PetscFunctionReturn(0);
2242d3f70b5SBarry Smith }
2252d3f70b5SBarry Smith 
2262d3f70b5SBarry Smith /*
2272d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2282d3f70b5SBarry Smith 
2292d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2302d3f70b5SBarry Smith */
2314a2ae208SSatish Balay #undef __FUNCT__
2324a2ae208SSatish Balay #define __FUNCT__ "TSPseudoJacobian"
2337bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2342d3f70b5SBarry Smith {
2352d3f70b5SBarry Smith   TS            ts = (TS) ctx;
2362d3f70b5SBarry Smith   int           ierr;
237ea709b57SSatish Balay   PetscScalar   mone = -1.0,mdt = 1.0/ts->time_step;
2382d3f70b5SBarry Smith 
2393a40ed3dSBarry Smith   PetscFunctionBegin;
2402d3f70b5SBarry Smith   /* construct users Jacobian */
241a7a1495cSBarry Smith   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
2422d3f70b5SBarry Smith 
243*614f654bSBarry Smith   /* shift and scale Jacobian */
2442d3f70b5SBarry Smith   ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
2452d3f70b5SBarry Smith   ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
246*614f654bSBarry Smith   if (*BB != *AA && *str != SAME_PRECONDITIONER) {
2472d3f70b5SBarry Smith     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
2482d3f70b5SBarry Smith     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
2492d3f70b5SBarry Smith   }
2502d3f70b5SBarry Smith 
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2522d3f70b5SBarry Smith }
2532d3f70b5SBarry Smith 
2542d3f70b5SBarry Smith 
2554a2ae208SSatish Balay #undef __FUNCT__
2564a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo"
2577bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts)
2582d3f70b5SBarry Smith {
2597bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
2602d3f70b5SBarry Smith   int       ierr,M,m;
2612d3f70b5SBarry Smith 
2623a40ed3dSBarry Smith   PetscFunctionBegin;
263b48991e4SBarry Smith   ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
2647bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2657bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2667bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2677bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
2692d3f70b5SBarry Smith }
2702d3f70b5SBarry Smith /*------------------------------------------------------------*/
2712d3f70b5SBarry Smith 
2724a2ae208SSatish Balay #undef __FUNCT__
2734a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultMonitor"
27487828ca2SBarry Smith int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
2752d3f70b5SBarry Smith {
2767bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
277d132466eSBarry Smith   int       ierr;
2782d3f70b5SBarry Smith 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
280142b95e3SSatish Balay   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
2822d3f70b5SBarry Smith }
2832d3f70b5SBarry Smith 
2844a2ae208SSatish Balay #undef __FUNCT__
2854a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo"
2867bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts)
2872d3f70b5SBarry Smith {
2884bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
289f1af5d2fSBarry Smith   int        ierr;
290f1af5d2fSBarry Smith   PetscTruth flg;
2912d3f70b5SBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2932d3f70b5SBarry Smith 
294b0a32e0cSBarry Smith   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
295b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
2962d3f70b5SBarry Smith     if (flg) {
297329f5518SBarry Smith       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
29828aa8177SBarry Smith     }
299b0a32e0cSBarry Smith     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
300ca4b7087SBarry Smith     if (flg) {
301ca4b7087SBarry Smith       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
302ca4b7087SBarry Smith     }
303419fbf4bSSatish Balay     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
304b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3052d3f70b5SBarry Smith 
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
3072d3f70b5SBarry Smith }
3082d3f70b5SBarry Smith 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo"
311b0a32e0cSBarry Smith static int TSView_Pseudo(TS ts,PetscViewer viewer)
3122d3f70b5SBarry Smith {
3133a40ed3dSBarry Smith   PetscFunctionBegin;
3143a40ed3dSBarry Smith   PetscFunctionReturn(0);
3152d3f70b5SBarry Smith }
3162d3f70b5SBarry Smith 
31782bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
3184a2ae208SSatish Balay #undef __FUNCT__
3194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
32082bf6240SBarry Smith /*@
32182bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
32282bf6240SBarry Smith    last timestep.
32382bf6240SBarry Smith 
32415091d37SBarry Smith    Collective on TS
32515091d37SBarry Smith 
32682bf6240SBarry Smith    Input Parameters:
32715091d37SBarry Smith +  ts - timestep context
32882bf6240SBarry Smith .  dt - user-defined function to verify timestep
32915091d37SBarry Smith -  ctx - [optional] user-defined context for private data
33082bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
33182bf6240SBarry Smith 
33215091d37SBarry Smith    Level: advanced
333fee21e36SBarry Smith 
33482bf6240SBarry Smith    Calling sequence of func:
33587828ca2SBarry Smith .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
33682bf6240SBarry Smith 
33782bf6240SBarry Smith .  update - latest solution vector
33882bf6240SBarry Smith .  ctx - [optional] timestep context
33982bf6240SBarry Smith .  newdt - the timestep to use for the next step
34082bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
34182bf6240SBarry Smith 
34282bf6240SBarry Smith    Notes:
34382bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
34482bf6240SBarry Smith    during the timestepping process.
34582bf6240SBarry Smith 
34682bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
34782bf6240SBarry Smith 
34882bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
34982bf6240SBarry Smith @*/
35087828ca2SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
35182bf6240SBarry Smith {
35287828ca2SBarry Smith   int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *);
35382bf6240SBarry Smith 
35482bf6240SBarry Smith   PetscFunctionBegin;
35582bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
35682bf6240SBarry Smith 
3571713a123SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)())&f);CHKERRQ(ierr);
35882bf6240SBarry Smith   if (f) {
35982bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
36082bf6240SBarry Smith   }
36182bf6240SBarry Smith   PetscFunctionReturn(0);
36282bf6240SBarry Smith }
36382bf6240SBarry Smith 
3644a2ae208SSatish Balay #undef __FUNCT__
3654a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
36682bf6240SBarry Smith /*@
36782bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
36882bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
36982bf6240SBarry Smith 
370fee21e36SBarry Smith    Collective on TS
371fee21e36SBarry Smith 
37215091d37SBarry Smith     Input Parameters:
37315091d37SBarry Smith +   ts - the timestep context
37415091d37SBarry Smith -   inc - the scaling factor >= 1.0
37515091d37SBarry Smith 
37682bf6240SBarry Smith     Options Database Key:
37782bf6240SBarry Smith $    -ts_pseudo_increment <increment>
37882bf6240SBarry Smith 
37915091d37SBarry Smith     Level: advanced
38015091d37SBarry Smith 
38182bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
38282bf6240SBarry Smith 
38382bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
38482bf6240SBarry Smith @*/
38587828ca2SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
38682bf6240SBarry Smith {
38787828ca2SBarry Smith   int ierr,(*f)(TS,PetscReal);
38882bf6240SBarry Smith 
38982bf6240SBarry Smith   PetscFunctionBegin;
39082bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
39182bf6240SBarry Smith 
3921713a123SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)())&f);CHKERRQ(ierr);
39382bf6240SBarry Smith   if (f) {
39482bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
39582bf6240SBarry Smith   }
39682bf6240SBarry Smith   PetscFunctionReturn(0);
39782bf6240SBarry Smith }
39882bf6240SBarry Smith 
3994a2ae208SSatish Balay #undef __FUNCT__
4004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
40182bf6240SBarry Smith /*@
40282bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
40382bf6240SBarry Smith     is computed via the formula
40482bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
40582bf6240SBarry Smith       rather than the default update,
40682bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
40782bf6240SBarry Smith 
40815091d37SBarry Smith    Collective on TS
40915091d37SBarry Smith 
41082bf6240SBarry Smith     Input Parameter:
41182bf6240SBarry Smith .   ts - the timestep context
41282bf6240SBarry Smith 
41382bf6240SBarry Smith     Options Database Key:
41482bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
41582bf6240SBarry Smith 
41615091d37SBarry Smith     Level: advanced
41715091d37SBarry Smith 
41882bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
41982bf6240SBarry Smith 
42082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
42182bf6240SBarry Smith @*/
42282bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts)
42382bf6240SBarry Smith {
42482bf6240SBarry Smith   int ierr,(*f)(TS);
42582bf6240SBarry Smith 
42682bf6240SBarry Smith   PetscFunctionBegin;
42782bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
42882bf6240SBarry Smith 
4291713a123SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)())&f);CHKERRQ(ierr);
43082bf6240SBarry Smith   if (f) {
43182bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
43282bf6240SBarry Smith   }
43382bf6240SBarry Smith   PetscFunctionReturn(0);
43482bf6240SBarry Smith }
43582bf6240SBarry Smith 
43682bf6240SBarry Smith 
4374a2ae208SSatish Balay #undef __FUNCT__
4384a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep"
43982bf6240SBarry Smith /*@
44082bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
44182bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
44282bf6240SBarry Smith 
44315091d37SBarry Smith    Collective on TS
44415091d37SBarry Smith 
44582bf6240SBarry Smith    Input Parameters:
44615091d37SBarry Smith +  ts - timestep context
44782bf6240SBarry Smith .  dt - function to compute timestep
44815091d37SBarry Smith -  ctx - [optional] user-defined context for private data
44982bf6240SBarry Smith          required by the function (may be PETSC_NULL)
45082bf6240SBarry Smith 
45115091d37SBarry Smith    Level: intermediate
452fee21e36SBarry Smith 
45382bf6240SBarry Smith    Calling sequence of func:
45487828ca2SBarry Smith .  func (TS ts,PetscReal *newdt,void *ctx);
45582bf6240SBarry Smith 
45682bf6240SBarry Smith .  newdt - the newly computed timestep
45782bf6240SBarry Smith .  ctx - [optional] timestep context
45882bf6240SBarry Smith 
45982bf6240SBarry Smith    Notes:
46082bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
46182bf6240SBarry Smith    during the timestepping process.
46282bf6240SBarry Smith 
46382bf6240SBarry Smith .keywords: timestep, pseudo, set
46482bf6240SBarry Smith 
46582bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
46682bf6240SBarry Smith @*/
46787828ca2SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
46882bf6240SBarry Smith {
46987828ca2SBarry Smith   int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *);
47082bf6240SBarry Smith 
47182bf6240SBarry Smith   PetscFunctionBegin;
47282bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
47382bf6240SBarry Smith 
4741713a123SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)())&f);CHKERRQ(ierr);
47582bf6240SBarry Smith   if (f) {
47682bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
47782bf6240SBarry Smith   }
47882bf6240SBarry Smith   PetscFunctionReturn(0);
47982bf6240SBarry Smith }
48082bf6240SBarry Smith 
48182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
48282bf6240SBarry Smith 
483fb2e594dSBarry Smith EXTERN_C_BEGIN
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
48687828ca2SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),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"
50187828ca2SBarry Smith int 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"
51482bf6240SBarry Smith int 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 
524fb2e594dSBarry Smith EXTERN_C_BEGIN
5254a2ae208SSatish Balay #undef __FUNCT__
5264a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
52787828ca2SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
52882bf6240SBarry Smith {
5294bbc92c1SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
53082bf6240SBarry Smith 
53182bf6240SBarry Smith   PetscFunctionBegin;
53282bf6240SBarry Smith   pseudo->dt      = dt;
53382bf6240SBarry Smith   pseudo->dtctx   = ctx;
53482bf6240SBarry Smith   PetscFunctionReturn(0);
53582bf6240SBarry Smith }
536fb2e594dSBarry Smith EXTERN_C_END
53782bf6240SBarry Smith 
53882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
53982bf6240SBarry Smith 
540fb2e594dSBarry Smith EXTERN_C_BEGIN
5414a2ae208SSatish Balay #undef __FUNCT__
5424a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo"
5437bf11e45SBarry Smith int TSCreate_Pseudo(TS ts)
5442d3f70b5SBarry Smith {
5457bf11e45SBarry Smith   TS_Pseudo  *pseudo;
5462d3f70b5SBarry Smith   int        ierr;
5472d3f70b5SBarry Smith 
5483a40ed3dSBarry Smith   PetscFunctionBegin;
5497bf11e45SBarry Smith   ts->destroy         = TSDestroy_Pseudo;
5507bf11e45SBarry Smith   ts->view            = TSView_Pseudo;
5512d3f70b5SBarry Smith 
5522d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
55329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
5542d3f70b5SBarry Smith   }
5552d3f70b5SBarry Smith   if (!ts->A) {
55629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
5572d3f70b5SBarry Smith   }
558273d9f13SBarry Smith 
5597bf11e45SBarry Smith   ts->setup           = TSSetUp_Pseudo;
5607bf11e45SBarry Smith   ts->step            = TSStep_Pseudo;
5617bf11e45SBarry Smith   ts->setfromoptions  = TSSetFromOptions_Pseudo;
5627bf11e45SBarry Smith 
5637bf11e45SBarry Smith   /* create the required nonlinear solver context */
5642d3f70b5SBarry Smith   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
5652d3f70b5SBarry Smith 
566b0a32e0cSBarry Smith   ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr);
567b0a32e0cSBarry Smith   PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
568eed86810SBarry Smith 
569549d3d68SSatish Balay   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
5707bf11e45SBarry Smith   ts->data = (void*)pseudo;
5712d3f70b5SBarry Smith 
57228aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
5734bbc92c1SBarry Smith   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
57428aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
57582bf6240SBarry Smith 
576f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
577e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
5780c97f09dSSatish Balay                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
579f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
580e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
5810c97f09dSSatish Balay                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
582f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
583e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
5840c97f09dSSatish Balay                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
5850c97f09dSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
5860c97f09dSSatish Balay                     "TSPseudoSetTimeStep_Pseudo",
5870c97f09dSSatish Balay                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
5883a40ed3dSBarry Smith   PetscFunctionReturn(0);
5892d3f70b5SBarry Smith }
590fb2e594dSBarry Smith EXTERN_C_END
5912d3f70b5SBarry Smith 
5924a2ae208SSatish Balay #undef __FUNCT__
5934a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep"
59482bf6240SBarry Smith /*@C
59582bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
59682bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
59728aa8177SBarry Smith 
59815091d37SBarry Smith    Collective on TS
59915091d37SBarry Smith 
60028aa8177SBarry Smith    Input Parameters:
60128aa8177SBarry Smith .  ts - the timestep context
60282bf6240SBarry Smith .  dtctx - unused timestep context
60328aa8177SBarry Smith 
60482bf6240SBarry Smith    Output Parameter:
60582bf6240SBarry Smith .  newdt - the timestep to use for the next step
60628aa8177SBarry Smith 
60715091d37SBarry Smith    Level: advanced
60815091d37SBarry Smith 
60982bf6240SBarry Smith .keywords: timestep, pseudo, default
610564e8f4eSLois Curfman McInnes 
61182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
61228aa8177SBarry Smith @*/
61387828ca2SBarry Smith int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
61428aa8177SBarry Smith {
61582bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
61687828ca2SBarry Smith   PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
61782bf6240SBarry Smith   int       ierr;
61828aa8177SBarry Smith 
6193a40ed3dSBarry Smith   PetscFunctionBegin;
62082bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
62182bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
62282bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
62382bf6240SBarry Smith     /* first time through so compute initial function norm */
62482bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
62582bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
62682bf6240SBarry Smith   }
62782bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
62882bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
629004884caSBarry Smith   } else if (pseudo->increment_dt_from_initial_dt) {
63082bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
63182bf6240SBarry Smith   } else {
63282bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
63382bf6240SBarry Smith   }
63482bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6353a40ed3dSBarry Smith   PetscFunctionReturn(0);
63628aa8177SBarry Smith }
6372d3f70b5SBarry Smith 
638004884caSBarry 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 
654