xref: /petsc/src/ts/impls/pseudo/posindep.c (revision f1af5d2ffeae1f5fc391a89939f4818e47770ae3)
1*f1af5d2fSBarry Smith /*$Id: posindep.c,v 1.37 1999/10/24 14:03:53 bsmith Exp bsmith $*/
22d3f70b5SBarry Smith /*
3fb4a63b6SLois Curfman McInnes        Code for Timestepping with implicit backwards Euler.
42d3f70b5SBarry Smith */
52d3f70b5SBarry Smith #include "src/ts/tsimpl.h"                /*I   "ts.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 
147bf11e45SBarry Smith   int    (*dt)(TS,double*,void*);              /* compute next timestep, and related context */
152d3f70b5SBarry Smith   void   *dtctx;
167bf11e45SBarry Smith   int    (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */
177bf11e45SBarry Smith   void   *verifyctx;
182d3f70b5SBarry Smith 
197bf11e45SBarry Smith   double initial_fnorm,fnorm;                  /* original and current norm of F(u) */
20ca4b7087SBarry Smith   double fnorm_previous;
2128aa8177SBarry Smith 
2228aa8177SBarry Smith   double dt_increment;                  /* scaling that dt is incremented each time-step */
23ca4b7087SBarry Smith   int    increment_dt_from_initial_dt;
247bf11e45SBarry Smith } TS_Pseudo;
252d3f70b5SBarry Smith 
262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/
272d3f70b5SBarry Smith 
2858562591SBarry Smith #undef __FUNC__
2958562591SBarry Smith #define __FUNC__ "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 @*/
527bf11e45SBarry Smith int TSPseudoComputeTimeStep(TS ts,double *dt)
537bf11e45SBarry Smith {
547bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
557bf11e45SBarry Smith   int       ierr;
567bf11e45SBarry Smith 
573a40ed3dSBarry Smith   PetscFunctionBegin;
587bf11e45SBarry Smith   PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
597bf11e45SBarry Smith   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
607bf11e45SBarry Smith   PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
613a40ed3dSBarry Smith   PetscFunctionReturn(0);
627bf11e45SBarry Smith }
637bf11e45SBarry Smith 
647bf11e45SBarry Smith 
657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/
6658562591SBarry Smith #undef __FUNC__
6758562591SBarry Smith #define __FUNC__ "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 @*/
927bf11e45SBarry Smith int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag)
937bf11e45SBarry Smith {
943a40ed3dSBarry Smith   PetscFunctionBegin;
957bf11e45SBarry Smith   *flag = 1;
963a40ed3dSBarry Smith   PetscFunctionReturn(0);
977bf11e45SBarry Smith }
987bf11e45SBarry Smith 
997bf11e45SBarry Smith 
10058562591SBarry Smith #undef __FUNC__
10158562591SBarry Smith #define __FUNC__ "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 @*/
1257bf11e45SBarry Smith int TSPseudoVerifyTimeStep(TS ts,Vec update,double *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 
14058562591SBarry Smith #undef __FUNC__
14158562591SBarry Smith #define __FUNC__ "TSStep_Pseudo"
1427bf11e45SBarry Smith static int TSStep_Pseudo(TS ts,int *steps,double *time)
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;
1477bf11e45SBarry Smith   double    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;
1722d3f70b5SBarry Smith   *time  = ts->ptime;
1733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1742d3f70b5SBarry Smith }
1752d3f70b5SBarry Smith 
1762d3f70b5SBarry Smith /*------------------------------------------------------------*/
17758562591SBarry Smith #undef __FUNC__
178d4bb536fSBarry Smith #define __FUNC__ "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);}
1882d3f70b5SBarry Smith   if (ts->Ashell)   {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
189606d414cSSatish Balay   ierr = PetscFree(pseudo);CHKERRQ(ierr);
1903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1912d3f70b5SBarry Smith }
1922d3f70b5SBarry Smith 
1932d3f70b5SBarry Smith 
1942d3f70b5SBarry Smith /*------------------------------------------------------------*/
1952d3f70b5SBarry Smith /*
1962d3f70b5SBarry Smith     This matrix shell multiply where user provided Shell matrix
1972d3f70b5SBarry Smith */
1982d3f70b5SBarry Smith 
19958562591SBarry Smith #undef __FUNC__
20058562591SBarry Smith #define __FUNC__ "TSPseudoMatMult"
2017bf11e45SBarry Smith int TSPseudoMatMult(Mat mat,Vec x,Vec y)
2022d3f70b5SBarry Smith {
2032d3f70b5SBarry Smith   TS     ts;
2042d3f70b5SBarry Smith   Scalar mdt,mone = -1.0;
2052d3f70b5SBarry Smith   int    ierr;
2062d3f70b5SBarry Smith 
2073a40ed3dSBarry Smith   PetscFunctionBegin;
2083a40ed3dSBarry Smith   ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr);
2092d3f70b5SBarry Smith   mdt = 1.0/ts->time_step;
2102d3f70b5SBarry Smith 
2112d3f70b5SBarry Smith   /* apply user provided function */
2122d3f70b5SBarry Smith   ierr = MatMult(ts->Ashell,x,y);CHKERRQ(ierr);
2132d3f70b5SBarry Smith   /* shift and scale by 1/dt - F */
2142d3f70b5SBarry Smith   ierr = VecAXPBY(&mdt,&mone,x,y);CHKERRQ(ierr);
2153a40ed3dSBarry Smith   PetscFunctionReturn(0);
2162d3f70b5SBarry Smith }
2172d3f70b5SBarry Smith 
2182d3f70b5SBarry Smith /*
2192d3f70b5SBarry Smith     This defines the nonlinear equation that is to be solved with SNES
2202d3f70b5SBarry Smith 
2212d3f70b5SBarry Smith               (U^{n+1} - U^{n})/dt - F(U^{n+1})
2222d3f70b5SBarry Smith */
22358562591SBarry Smith #undef __FUNC__
22458562591SBarry Smith #define __FUNC__ "TSPseudoFunction"
2257bf11e45SBarry Smith int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
2262d3f70b5SBarry Smith {
2272d3f70b5SBarry Smith   TS     ts = (TS) ctx;
2282d3f70b5SBarry Smith   Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
2292d3f70b5SBarry Smith   int    ierr,i,n;
2302d3f70b5SBarry Smith 
2313a40ed3dSBarry Smith   PetscFunctionBegin;
2322d3f70b5SBarry Smith   /* apply user provided function */
2332d3f70b5SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
2347bf11e45SBarry Smith   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
2352d3f70b5SBarry Smith   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
2362d3f70b5SBarry Smith   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
2372d3f70b5SBarry Smith   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
2382d3f70b5SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
2392d3f70b5SBarry Smith   for ( i=0; i<n; i++ ) {
2402d3f70b5SBarry Smith     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
2412d3f70b5SBarry Smith   }
242888f2ed8SSatish Balay   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
243888f2ed8SSatish Balay   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
244888f2ed8SSatish Balay   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
2452d3f70b5SBarry Smith 
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
2472d3f70b5SBarry Smith }
2482d3f70b5SBarry Smith 
2492d3f70b5SBarry Smith /*
2502d3f70b5SBarry Smith    This constructs the Jacobian needed for SNES
2512d3f70b5SBarry Smith 
2522d3f70b5SBarry Smith              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
2532d3f70b5SBarry Smith */
25458562591SBarry Smith #undef __FUNC__
25558562591SBarry Smith #define __FUNC__ "TSPseudoJacobian"
2567bf11e45SBarry Smith int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
2572d3f70b5SBarry Smith {
2582d3f70b5SBarry Smith   TS      ts = (TS) ctx;
2592d3f70b5SBarry Smith   int     ierr;
2602d3f70b5SBarry Smith   Scalar  mone = -1.0, mdt = 1.0/ts->time_step;
2612d3f70b5SBarry Smith   MatType mtype;
2622d3f70b5SBarry Smith 
2633a40ed3dSBarry Smith   PetscFunctionBegin;
2642d3f70b5SBarry Smith   /* construct users Jacobian */
2652d3f70b5SBarry Smith   if (ts->rhsjacobian) {
2662d3f70b5SBarry Smith     ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr);
2672d3f70b5SBarry Smith   }
2682d3f70b5SBarry Smith 
2692d3f70b5SBarry Smith   /* shift and scale Jacobian, if not a shell matrix */
270888f2ed8SSatish Balay   ierr = MatGetType(*AA,&mtype,PETSC_NULL);CHKERRQ(ierr);
2712d3f70b5SBarry Smith   if (mtype != MATSHELL) {
2722d3f70b5SBarry Smith     ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
2732d3f70b5SBarry Smith     ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
2742d3f70b5SBarry Smith   }
275888f2ed8SSatish Balay   ierr = MatGetType(*BB,&mtype,PETSC_NULL);CHKERRQ(ierr);
2762d3f70b5SBarry Smith   if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) {
2772d3f70b5SBarry Smith     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
2782d3f70b5SBarry Smith     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
2792d3f70b5SBarry Smith   }
2802d3f70b5SBarry Smith 
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
2822d3f70b5SBarry Smith }
2832d3f70b5SBarry Smith 
2842d3f70b5SBarry Smith 
28558562591SBarry Smith #undef __FUNC__
28658562591SBarry Smith #define __FUNC__ "TSSetUp_Pseudo"
2877bf11e45SBarry Smith static int TSSetUp_Pseudo(TS ts)
2882d3f70b5SBarry Smith {
2897bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
2902d3f70b5SBarry Smith   int       ierr, M, m;
2912d3f70b5SBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2937bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
2947bf11e45SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
2957bf11e45SBarry Smith   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
2962d3f70b5SBarry Smith   if (ts->Ashell) { /* construct new shell matrix */
2972d3f70b5SBarry Smith     ierr = VecGetSize(ts->vec_sol,&M);CHKERRQ(ierr);
2982d3f70b5SBarry Smith     ierr = VecGetLocalSize(ts->vec_sol,&m);CHKERRQ(ierr);
2992d3f70b5SBarry Smith     ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A);CHKERRQ(ierr);
3007bf11e45SBarry Smith     ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr);
3012d3f70b5SBarry Smith   }
3027bf11e45SBarry Smith   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
3033a40ed3dSBarry Smith   PetscFunctionReturn(0);
3042d3f70b5SBarry Smith }
3052d3f70b5SBarry Smith /*------------------------------------------------------------*/
3062d3f70b5SBarry Smith 
30758562591SBarry Smith #undef __FUNC__
308d4bb536fSBarry Smith #define __FUNC__ "TSPseudoDefaultMonitor"
3092d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
3102d3f70b5SBarry Smith {
3117bf11e45SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
312d132466eSBarry Smith   int       ierr;
3132d3f70b5SBarry Smith 
3143a40ed3dSBarry Smith   PetscFunctionBegin;
315d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm);CHKERRQ(ierr);
3163a40ed3dSBarry Smith   PetscFunctionReturn(0);
3172d3f70b5SBarry Smith }
3182d3f70b5SBarry Smith 
31958562591SBarry Smith #undef __FUNC__
32058562591SBarry Smith #define __FUNC__ "TSSetFromOptions_Pseudo"
3217bf11e45SBarry Smith static int TSSetFromOptions_Pseudo(TS ts)
3222d3f70b5SBarry Smith {
323*f1af5d2fSBarry Smith   int        ierr;
324*f1af5d2fSBarry Smith   PetscTruth flg;
32528aa8177SBarry Smith   double     inc;
3262d3f70b5SBarry Smith 
3273a40ed3dSBarry Smith   PetscFunctionBegin;
3282d3f70b5SBarry Smith   ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
3292d3f70b5SBarry Smith 
3302d3f70b5SBarry Smith   ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg);CHKERRQ(ierr);
3312d3f70b5SBarry Smith   if (flg) {
33228aa8177SBarry Smith     ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0);CHKERRQ(ierr);
33328aa8177SBarry Smith   }
33428aa8177SBarry Smith   ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg);CHKERRQ(ierr);
33528aa8177SBarry Smith   if (flg) {
33628aa8177SBarry Smith     ierr = TSPseudoSetTimeStepIncrement(ts,inc);CHKERRQ(ierr);
3372d3f70b5SBarry Smith   }
338ca4b7087SBarry Smith   ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr);
339ca4b7087SBarry Smith   if (flg) {
340ca4b7087SBarry Smith     ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
341ca4b7087SBarry Smith   }
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
3432d3f70b5SBarry Smith }
3442d3f70b5SBarry Smith 
34558562591SBarry Smith #undef __FUNC__
346d4bb536fSBarry Smith #define __FUNC__ "TSPrintHelp_Pseudo"
347ca4b7087SBarry Smith static int TSPrintHelp_Pseudo(TS ts,char *p)
3482d3f70b5SBarry Smith {
349d132466eSBarry Smith   int ierr;
350d132466eSBarry Smith 
3513a40ed3dSBarry Smith   PetscFunctionBegin;
352d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n");CHKERRQ(ierr);
353d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p);CHKERRQ(ierr);
354d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p);CHKERRQ(ierr);
355d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm,"     initial fnorm/current fnorm to determine new timestep\n");CHKERRQ(ierr);
356d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(ts->comm,"     default is current_dt * previous fnorm/current fnorm\n");CHKERRQ(ierr);
3573a40ed3dSBarry Smith   PetscFunctionReturn(0);
3582d3f70b5SBarry Smith }
3592d3f70b5SBarry Smith 
36058562591SBarry Smith #undef __FUNC__
361d4bb536fSBarry Smith #define __FUNC__ "TSView_Pseudo"
362e1311b90SBarry Smith static int TSView_Pseudo(TS ts,Viewer viewer)
3632d3f70b5SBarry Smith {
3643a40ed3dSBarry Smith   PetscFunctionBegin;
3653a40ed3dSBarry Smith   PetscFunctionReturn(0);
3662d3f70b5SBarry Smith }
3672d3f70b5SBarry Smith 
36882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
36982bf6240SBarry Smith #undef __FUNC__
37082bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep"
37182bf6240SBarry Smith /*@
37282bf6240SBarry Smith    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
37382bf6240SBarry Smith    last timestep.
37482bf6240SBarry Smith 
37515091d37SBarry Smith    Collective on TS
37615091d37SBarry Smith 
37782bf6240SBarry Smith    Input Parameters:
37815091d37SBarry Smith +  ts - timestep context
37982bf6240SBarry Smith .  dt - user-defined function to verify timestep
38015091d37SBarry Smith -  ctx - [optional] user-defined context for private data
38182bf6240SBarry Smith          for the timestep verification routine (may be PETSC_NULL)
38282bf6240SBarry Smith 
38315091d37SBarry Smith    Level: advanced
384fee21e36SBarry Smith 
38582bf6240SBarry Smith    Calling sequence of func:
38682bf6240SBarry Smith .  func (TS ts,Vec update,void *ctx,double *newdt,int *flag);
38782bf6240SBarry Smith 
38882bf6240SBarry Smith .  update - latest solution vector
38982bf6240SBarry Smith .  ctx - [optional] timestep context
39082bf6240SBarry Smith .  newdt - the timestep to use for the next step
39182bf6240SBarry Smith .  flag - flag indicating whether the last time step was acceptable
39282bf6240SBarry Smith 
39382bf6240SBarry Smith    Notes:
39482bf6240SBarry Smith    The routine set here will be called by TSPseudoVerifyTimeStep()
39582bf6240SBarry Smith    during the timestepping process.
39682bf6240SBarry Smith 
39782bf6240SBarry Smith .keywords: timestep, pseudo, set, verify
39882bf6240SBarry Smith 
39982bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
40082bf6240SBarry Smith @*/
40182bf6240SBarry Smith int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
40282bf6240SBarry Smith {
40382bf6240SBarry Smith   int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *);
40482bf6240SBarry Smith 
40582bf6240SBarry Smith   PetscFunctionBegin;
40682bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
40782bf6240SBarry Smith 
40894d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void **)&f);CHKERRQ(ierr);
40982bf6240SBarry Smith   if (f) {
41082bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
41182bf6240SBarry Smith   }
41282bf6240SBarry Smith   PetscFunctionReturn(0);
41382bf6240SBarry Smith }
41482bf6240SBarry Smith 
41582bf6240SBarry Smith #undef __FUNC__
41682bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement"
41782bf6240SBarry Smith /*@
41882bf6240SBarry Smith     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
41982bf6240SBarry Smith     dt when using the TSPseudoDefaultTimeStep() routine.
42082bf6240SBarry Smith 
421fee21e36SBarry Smith    Collective on TS
422fee21e36SBarry Smith 
42315091d37SBarry Smith     Input Parameters:
42415091d37SBarry Smith +   ts - the timestep context
42515091d37SBarry Smith -   inc - the scaling factor >= 1.0
42615091d37SBarry Smith 
42782bf6240SBarry Smith     Options Database Key:
42882bf6240SBarry Smith $    -ts_pseudo_increment <increment>
42982bf6240SBarry Smith 
43015091d37SBarry Smith     Level: advanced
43115091d37SBarry Smith 
43282bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
43382bf6240SBarry Smith 
43482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
43582bf6240SBarry Smith @*/
43682bf6240SBarry Smith int TSPseudoSetTimeStepIncrement(TS ts,double inc)
43782bf6240SBarry Smith {
43882bf6240SBarry Smith   int ierr, (*f)(TS,double);
43982bf6240SBarry Smith 
44082bf6240SBarry Smith   PetscFunctionBegin;
44182bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
44282bf6240SBarry Smith 
44394d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void **)&f);CHKERRQ(ierr);
44482bf6240SBarry Smith   if (f) {
44582bf6240SBarry Smith     ierr = (*f)(ts,inc);CHKERRQ(ierr);
44682bf6240SBarry Smith   }
44782bf6240SBarry Smith   PetscFunctionReturn(0);
44882bf6240SBarry Smith }
44982bf6240SBarry Smith 
45082bf6240SBarry Smith #undef __FUNC__
45182bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt"
45282bf6240SBarry Smith /*@
45382bf6240SBarry Smith     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
45482bf6240SBarry Smith     is computed via the formula
45582bf6240SBarry Smith $         dt = initial_dt*initial_fnorm/current_fnorm
45682bf6240SBarry Smith       rather than the default update,
45782bf6240SBarry Smith $         dt = current_dt*previous_fnorm/current_fnorm.
45882bf6240SBarry Smith 
45915091d37SBarry Smith    Collective on TS
46015091d37SBarry Smith 
46182bf6240SBarry Smith     Input Parameter:
46282bf6240SBarry Smith .   ts - the timestep context
46382bf6240SBarry Smith 
46482bf6240SBarry Smith     Options Database Key:
46582bf6240SBarry Smith $    -ts_pseudo_increment_dt_from_initial_dt
46682bf6240SBarry Smith 
46715091d37SBarry Smith     Level: advanced
46815091d37SBarry Smith 
46982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment
47082bf6240SBarry Smith 
47182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
47282bf6240SBarry Smith @*/
47382bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt(TS ts)
47482bf6240SBarry Smith {
47582bf6240SBarry Smith   int ierr, (*f)(TS);
47682bf6240SBarry Smith 
47782bf6240SBarry Smith   PetscFunctionBegin;
47882bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
47982bf6240SBarry Smith 
48094d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void **)&f);CHKERRQ(ierr);
48182bf6240SBarry Smith   if (f) {
48282bf6240SBarry Smith     ierr = (*f)(ts);CHKERRQ(ierr);
48382bf6240SBarry Smith   }
48482bf6240SBarry Smith   PetscFunctionReturn(0);
48582bf6240SBarry Smith }
48682bf6240SBarry Smith 
48782bf6240SBarry Smith 
48882bf6240SBarry Smith #undef __FUNC__
48982bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep"
49082bf6240SBarry Smith /*@
49182bf6240SBarry Smith    TSPseudoSetTimeStep - Sets the user-defined routine to be
49282bf6240SBarry Smith    called at each pseudo-timestep to update the timestep.
49382bf6240SBarry Smith 
49415091d37SBarry Smith    Collective on TS
49515091d37SBarry Smith 
49682bf6240SBarry Smith    Input Parameters:
49715091d37SBarry Smith +  ts - timestep context
49882bf6240SBarry Smith .  dt - function to compute timestep
49915091d37SBarry Smith -  ctx - [optional] user-defined context for private data
50082bf6240SBarry Smith          required by the function (may be PETSC_NULL)
50182bf6240SBarry Smith 
50215091d37SBarry Smith    Level: intermediate
503fee21e36SBarry Smith 
50482bf6240SBarry Smith    Calling sequence of func:
50582bf6240SBarry Smith .  func (TS ts,double *newdt,void *ctx);
50682bf6240SBarry Smith 
50782bf6240SBarry Smith .  newdt - the newly computed timestep
50882bf6240SBarry Smith .  ctx - [optional] timestep context
50982bf6240SBarry Smith 
51082bf6240SBarry Smith    Notes:
51182bf6240SBarry Smith    The routine set here will be called by TSPseudoComputeTimeStep()
51282bf6240SBarry Smith    during the timestepping process.
51382bf6240SBarry Smith 
51482bf6240SBarry Smith .keywords: timestep, pseudo, set
51582bf6240SBarry Smith 
51682bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
51782bf6240SBarry Smith @*/
51882bf6240SBarry Smith int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx)
51982bf6240SBarry Smith {
52082bf6240SBarry Smith   int ierr, (*f)(TS,int (*)(TS,double *,void *),void *);
52182bf6240SBarry Smith 
52282bf6240SBarry Smith   PetscFunctionBegin;
52382bf6240SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
52482bf6240SBarry Smith 
52594d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void **)&f);CHKERRQ(ierr);
52682bf6240SBarry Smith   if (f) {
52782bf6240SBarry Smith     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
52882bf6240SBarry Smith   }
52982bf6240SBarry Smith   PetscFunctionReturn(0);
53082bf6240SBarry Smith }
53182bf6240SBarry Smith 
53282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
53382bf6240SBarry Smith 
534fb2e594dSBarry Smith EXTERN_C_BEGIN
53582bf6240SBarry Smith #undef __FUNC__
53682bf6240SBarry Smith #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo"
53782bf6240SBarry Smith int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
53882bf6240SBarry Smith {
53982bf6240SBarry Smith   TS_Pseudo *pseudo;
54082bf6240SBarry Smith 
54182bf6240SBarry Smith   PetscFunctionBegin;
54282bf6240SBarry Smith   pseudo              = (TS_Pseudo*) ts->data;
54382bf6240SBarry Smith   pseudo->verify      = dt;
54482bf6240SBarry Smith   pseudo->verifyctx   = ctx;
54582bf6240SBarry Smith   PetscFunctionReturn(0);
54682bf6240SBarry Smith }
547fb2e594dSBarry Smith EXTERN_C_END
54882bf6240SBarry Smith 
549fb2e594dSBarry Smith EXTERN_C_BEGIN
55082bf6240SBarry Smith #undef __FUNC__
55182bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo"
55282bf6240SBarry Smith int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc)
55382bf6240SBarry Smith {
55482bf6240SBarry Smith   TS_Pseudo *pseudo;
55582bf6240SBarry Smith 
55682bf6240SBarry Smith   PetscFunctionBegin;
55782bf6240SBarry Smith   pseudo               = (TS_Pseudo*) ts->data;
55882bf6240SBarry Smith   pseudo->dt_increment = inc;
55982bf6240SBarry Smith   PetscFunctionReturn(0);
56082bf6240SBarry Smith }
561fb2e594dSBarry Smith EXTERN_C_END
56282bf6240SBarry Smith 
563fb2e594dSBarry Smith EXTERN_C_BEGIN
56482bf6240SBarry Smith #undef __FUNC__
56582bf6240SBarry Smith #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
56682bf6240SBarry Smith int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
56782bf6240SBarry Smith {
56882bf6240SBarry Smith   TS_Pseudo *pseudo;
56982bf6240SBarry Smith 
57082bf6240SBarry Smith   PetscFunctionBegin;
57182bf6240SBarry Smith   pseudo                               = (TS_Pseudo*) ts->data;
57282bf6240SBarry Smith   pseudo->increment_dt_from_initial_dt = 1;
57382bf6240SBarry Smith   PetscFunctionReturn(0);
57482bf6240SBarry Smith }
575fb2e594dSBarry Smith EXTERN_C_END
57682bf6240SBarry Smith 
577fb2e594dSBarry Smith EXTERN_C_BEGIN
57882bf6240SBarry Smith #undef __FUNC__
57982bf6240SBarry Smith #define __FUNC__ "TSPseudoSetTimeStep_Pseudo"
58082bf6240SBarry Smith int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx)
58182bf6240SBarry Smith {
58282bf6240SBarry Smith   TS_Pseudo *pseudo;
58382bf6240SBarry Smith 
58482bf6240SBarry Smith   PetscFunctionBegin;
58582bf6240SBarry Smith   pseudo          = (TS_Pseudo*) ts->data;
58682bf6240SBarry Smith   pseudo->dt      = dt;
58782bf6240SBarry Smith   pseudo->dtctx   = ctx;
58882bf6240SBarry Smith   PetscFunctionReturn(0);
58982bf6240SBarry Smith }
590fb2e594dSBarry Smith EXTERN_C_END
59182bf6240SBarry Smith 
59282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */
59382bf6240SBarry Smith 
594fb2e594dSBarry Smith EXTERN_C_BEGIN
59558562591SBarry Smith #undef __FUNC__
59658562591SBarry Smith #define __FUNC__ "TSCreate_Pseudo"
5977bf11e45SBarry Smith int TSCreate_Pseudo(TS ts )
5982d3f70b5SBarry Smith {
5997bf11e45SBarry Smith   TS_Pseudo *pseudo;
6002d3f70b5SBarry Smith   int       ierr;
6012d3f70b5SBarry Smith   MatType   mtype;
6022d3f70b5SBarry Smith 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
6047bf11e45SBarry Smith   ts->destroy         = TSDestroy_Pseudo;
6057bf11e45SBarry Smith   ts->printhelp       = TSPrintHelp_Pseudo;
6067bf11e45SBarry Smith   ts->view            = TSView_Pseudo;
6072d3f70b5SBarry Smith 
6082d3f70b5SBarry Smith   if (ts->problem_type == TS_LINEAR) {
609a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems");
6102d3f70b5SBarry Smith   }
6112d3f70b5SBarry Smith   if (!ts->A) {
612a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian");
6132d3f70b5SBarry Smith   }
614888f2ed8SSatish Balay   ierr = MatGetType(ts->A,&mtype,PETSC_NULL);CHKERRQ(ierr);
6152d3f70b5SBarry Smith   if (mtype == MATSHELL) {
6162d3f70b5SBarry Smith     ts->Ashell = ts->A;
6172d3f70b5SBarry Smith   }
6187bf11e45SBarry Smith   ts->setup           = TSSetUp_Pseudo;
6197bf11e45SBarry Smith   ts->step            = TSStep_Pseudo;
6207bf11e45SBarry Smith   ts->setfromoptions  = TSSetFromOptions_Pseudo;
6217bf11e45SBarry Smith 
6227bf11e45SBarry Smith   /* create the required nonlinear solver context */
6232d3f70b5SBarry Smith   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
6242d3f70b5SBarry Smith 
6257bf11e45SBarry Smith   pseudo   = PetscNew(TS_Pseudo);CHKPTRQ(pseudo);
626eed86810SBarry Smith   PLogObjectMemory(ts,sizeof(TS_Pseudo));
627eed86810SBarry Smith 
628549d3d68SSatish Balay   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
6297bf11e45SBarry Smith   ts->data = (void *) pseudo;
6302d3f70b5SBarry Smith 
63128aa8177SBarry Smith   pseudo->dt_increment                 = 1.1;
632ca4b7087SBarry Smith   pseudo->increment_dt_from_initial_dt = 0;
63328aa8177SBarry Smith   pseudo->dt                           = TSPseudoDefaultTimeStep;
63482bf6240SBarry Smith 
635*f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
636e1311b90SBarry Smith                     "TSPseudoSetVerifyTimeStep_Pseudo",
637e1311b90SBarry Smith                     (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
638*f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
639e1311b90SBarry Smith                     "TSPseudoSetTimeStepIncrement_Pseudo",
640e1311b90SBarry Smith                     (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
641*f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
642e1311b90SBarry Smith                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
643e1311b90SBarry Smith                     (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
644*f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","TSPseudoSetTimeStep_Pseudo",
645e1311b90SBarry Smith                     (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
6463a40ed3dSBarry Smith   PetscFunctionReturn(0);
6472d3f70b5SBarry Smith }
648fb2e594dSBarry Smith EXTERN_C_END
6492d3f70b5SBarry Smith 
65058562591SBarry Smith #undef __FUNC__
65182bf6240SBarry Smith #define __FUNC__ "TSPseudoDefaultTimeStep"
65282bf6240SBarry Smith /*@C
65382bf6240SBarry Smith    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
65482bf6240SBarry Smith    Use with TSPseudoSetTimeStep().
65528aa8177SBarry Smith 
65615091d37SBarry Smith    Collective on TS
65715091d37SBarry Smith 
65828aa8177SBarry Smith    Input Parameters:
65928aa8177SBarry Smith .  ts - the timestep context
66082bf6240SBarry Smith .  dtctx - unused timestep context
66128aa8177SBarry Smith 
66282bf6240SBarry Smith    Output Parameter:
66382bf6240SBarry Smith .  newdt - the timestep to use for the next step
66428aa8177SBarry Smith 
66515091d37SBarry Smith    Level: advanced
66615091d37SBarry Smith 
66782bf6240SBarry Smith .keywords: timestep, pseudo, default
668564e8f4eSLois Curfman McInnes 
66982bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
67028aa8177SBarry Smith @*/
67182bf6240SBarry Smith int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx)
67228aa8177SBarry Smith {
67382bf6240SBarry Smith   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
67482bf6240SBarry Smith   double    inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
67582bf6240SBarry Smith   int       ierr;
67628aa8177SBarry Smith 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
67882bf6240SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
67982bf6240SBarry Smith   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
68082bf6240SBarry Smith   if (pseudo->initial_fnorm == 0.0) {
68182bf6240SBarry Smith     /* first time through so compute initial function norm */
68282bf6240SBarry Smith     pseudo->initial_fnorm = pseudo->fnorm;
68382bf6240SBarry Smith     fnorm_previous        = pseudo->fnorm;
68482bf6240SBarry Smith   }
68582bf6240SBarry Smith   if (pseudo->fnorm == 0.0) {
68682bf6240SBarry Smith     *newdt = 1.e12*inc*ts->time_step;
68782bf6240SBarry Smith   }
68882bf6240SBarry Smith   else if (pseudo->increment_dt_from_initial_dt) {
68982bf6240SBarry Smith     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
69082bf6240SBarry Smith   } else {
69182bf6240SBarry Smith     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
69282bf6240SBarry Smith   }
69382bf6240SBarry Smith   pseudo->fnorm_previous = pseudo->fnorm;
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
69528aa8177SBarry Smith }
6962d3f70b5SBarry Smith 
697