xref: /petsc/src/ts/impls/pseudo/posindep.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 /*
2        Code for Timestepping with implicit backwards Euler.
3 */
4 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5 
6 typedef struct {
7   Vec update;       /* work vector where new solution is formed */
8   Vec func;         /* work vector where F(t[i],u[i]) is stored */
9   Vec xdot;         /* work vector for time derivative of state */
10 
11   /* information used for Pseudo-timestepping */
12 
13   PetscErrorCode (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
14   void *dtctx;
15   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscBool*);  /* verify previous timestep and related context */
16   void *verifyctx;
17 
18   PetscReal fnorm_initial,fnorm;                   /* original and current norm of F(u) */
19   PetscReal fnorm_previous;
20 
21   PetscReal dt_initial;                     /* initial time-step */
22   PetscReal dt_increment;                   /* scaling that dt is incremented each time-step */
23   PetscReal dt_max;                         /* maximum time step */
24   PetscBool increment_dt_from_initial_dt;
25   PetscReal fatol,frtol;
26 } TS_Pseudo;
27 
28 /* ------------------------------------------------------------------------------*/
29 
30 /*@C
31     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32     pseudo-timestepping process.
33 
34     Collective on TS
35 
36     Input Parameter:
37 .   ts - timestep context
38 
39     Output Parameter:
40 .   dt - newly computed timestep
41 
42     Level: developer
43 
44     Notes:
45     The routine to be called here to compute the timestep should be
46     set by calling TSPseudoSetTimeStep().
47 
48 .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()`
49 @*/
50 PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
51 {
52   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
53 
54   PetscFunctionBegin;
55   PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0));
56   PetscCall((*pseudo->dt)(ts,dt,pseudo->dtctx));
57   PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0));
58   PetscFunctionReturn(0);
59 }
60 
61 /* ------------------------------------------------------------------------------*/
62 /*@C
63    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
64 
65    Collective on TS
66 
67    Input Parameters:
68 +  ts - the timestep context
69 .  dtctx - unused timestep context
70 -  update - latest solution vector
71 
72    Output Parameters:
73 +  newdt - the timestep to use for the next step
74 -  flag - flag indicating whether the last time step was acceptable
75 
76    Level: advanced
77 
78    Note:
79    This routine always returns a flag of 1, indicating an acceptable
80    timestep.
81 
82 .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()`
83 @*/
84 PetscErrorCode  TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
85 {
86   PetscFunctionBegin;
87   *flag = PETSC_TRUE;
88   PetscFunctionReturn(0);
89 }
90 
91 /*@
92     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
93 
94     Collective on TS
95 
96     Input Parameters:
97 +   ts - timestep context
98 -   update - latest solution vector
99 
100     Output Parameters:
101 +   dt - newly computed timestep (if it had to shrink)
102 -   flag - indicates if current timestep was ok
103 
104     Level: advanced
105 
106     Notes:
107     The routine to be called here to compute the timestep should be
108     set by calling TSPseudoSetVerifyTimeStep().
109 
110 .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()`
111 @*/
112 PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag)
113 {
114   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
115 
116   PetscFunctionBegin;
117   *flag = PETSC_TRUE;
118   if (pseudo->verify) PetscCall((*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag));
119   PetscFunctionReturn(0);
120 }
121 
122 /* --------------------------------------------------------------------------------*/
123 
124 static PetscErrorCode TSStep_Pseudo(TS ts)
125 {
126   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
127   PetscInt            nits,lits,reject;
128   PetscBool           stepok;
129   PetscReal           next_time_step = ts->time_step;
130 
131   PetscFunctionBegin;
132   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
133   PetscCall(VecCopy(ts->vec_sol,pseudo->update));
134   PetscCall(TSPseudoComputeTimeStep(ts,&next_time_step));
135   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
136     ts->time_step = next_time_step;
137     PetscCall(TSPreStage(ts,ts->ptime+ts->time_step));
138     PetscCall(SNESSolve(ts->snes,NULL,pseudo->update));
139     PetscCall(SNESGetIterationNumber(ts->snes,&nits));
140     PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits));
141     ts->snes_its += nits; ts->ksp_its += lits;
142     PetscCall(TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update)));
143     PetscCall(TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok));
144     if (!stepok) {next_time_step = ts->time_step; continue;}
145     pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
146     PetscCall(TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok));
147     if (stepok) break;
148   }
149   if (reject >= ts->max_reject) {
150     ts->reason = TS_DIVERGED_STEP_REJECTED;
151     PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,reject));
152     PetscFunctionReturn(0);
153   }
154 
155   PetscCall(VecCopy(pseudo->update,ts->vec_sol));
156   ts->ptime += ts->time_step;
157   ts->time_step = next_time_step;
158 
159   if (pseudo->fnorm < 0) {
160     PetscCall(VecZeroEntries(pseudo->xdot));
161     PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
162     PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
163   }
164   if (pseudo->fnorm < pseudo->fatol) {
165     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
166     PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n",ts->steps,(double)pseudo->fnorm,(double)pseudo->frtol));
167     PetscFunctionReturn(0);
168   }
169   if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) {
170     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
171     PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,(double)pseudo->fnorm,(double)pseudo->fnorm_initial,(double)pseudo->fatol));
172     PetscFunctionReturn(0);
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 /*------------------------------------------------------------*/
178 static PetscErrorCode TSReset_Pseudo(TS ts)
179 {
180   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
181 
182   PetscFunctionBegin;
183   PetscCall(VecDestroy(&pseudo->update));
184   PetscCall(VecDestroy(&pseudo->func));
185   PetscCall(VecDestroy(&pseudo->xdot));
186   PetscFunctionReturn(0);
187 }
188 
189 static PetscErrorCode TSDestroy_Pseudo(TS ts)
190 {
191   PetscFunctionBegin;
192   PetscCall(TSReset_Pseudo(ts));
193   PetscCall(PetscFree(ts->data));
194   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL));
195   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL));
196   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL));
197   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL));
198   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL));
199   PetscFunctionReturn(0);
200 }
201 
202 /*------------------------------------------------------------*/
203 
204 /*
205     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
206 */
207 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
208 {
209   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
210   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
211   PetscScalar       *xdot;
212   PetscInt          i,n;
213 
214   PetscFunctionBegin;
215   *Xdot = NULL;
216   PetscCall(VecGetArrayRead(ts->vec_sol,&xn));
217   PetscCall(VecGetArrayRead(X,&xnp1));
218   PetscCall(VecGetArray(pseudo->xdot,&xdot));
219   PetscCall(VecGetLocalSize(X,&n));
220   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
221   PetscCall(VecRestoreArrayRead(ts->vec_sol,&xn));
222   PetscCall(VecRestoreArrayRead(X,&xnp1));
223   PetscCall(VecRestoreArray(pseudo->xdot,&xdot));
224   *Xdot = pseudo->xdot;
225   PetscFunctionReturn(0);
226 }
227 
228 /*
229     The transient residual is
230 
231         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
232 
233     or for ODE,
234 
235         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
236 
237     This is the function that must be evaluated for transient simulation and for
238     finite difference Jacobians.  On the first Newton step, this algorithm uses
239     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
240     residual is actually the steady state residual.  Pseudotransient
241     continuation as described in the literature is a linearly implicit
242     algorithm, it only takes this one Newton step with the steady state
243     residual, and then advances to the next time step.
244 */
245 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
246 {
247   Vec            Xdot;
248 
249   PetscFunctionBegin;
250   PetscCall(TSPseudoGetXdot(ts,X,&Xdot));
251   PetscCall(TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE));
252   PetscFunctionReturn(0);
253 }
254 
255 /*
256    This constructs the Jacobian needed for SNES.  For DAE, this is
257 
258        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
259 
260     and for ODE:
261 
262        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
263 */
264 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
265 {
266   Vec            Xdot;
267 
268   PetscFunctionBegin;
269   PetscCall(TSPseudoGetXdot(ts,X,&Xdot));
270   PetscCall(TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE));
271   PetscFunctionReturn(0);
272 }
273 
274 static PetscErrorCode TSSetUp_Pseudo(TS ts)
275 {
276   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
277 
278   PetscFunctionBegin;
279   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->update));
280   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->func));
281   PetscCall(VecDuplicate(ts->vec_sol,&pseudo->xdot));
282   PetscFunctionReturn(0);
283 }
284 /*------------------------------------------------------------*/
285 
286 static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
287 {
288   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
289   PetscViewer    viewer = (PetscViewer) dummy;
290 
291   PetscFunctionBegin;
292   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
293     PetscCall(VecZeroEntries(pseudo->xdot));
294     PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
295     PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
296   }
297   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel));
298   PetscCall(PetscViewerASCIIPrintf(viewer,"TS %" PetscInt_FMT " dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm));
299   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel));
300   PetscFunctionReturn(0);
301 }
302 
303 static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts)
304 {
305   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
306   PetscBool      flg = PETSC_FALSE;
307   PetscViewer    viewer;
308 
309   PetscFunctionBegin;
310   PetscOptionsHeadBegin(PetscOptionsObject,"Pseudo-timestepping options");
311   PetscCall(PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL));
312   if (flg) {
313     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer));
314     PetscCall(TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
315   }
316   flg  = pseudo->increment_dt_from_initial_dt;
317   PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL));
318   pseudo->increment_dt_from_initial_dt = flg;
319   PetscCall(PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL));
320   PetscCall(PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL));
321   PetscCall(PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL));
322   PetscCall(PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL));
323   PetscOptionsHeadEnd();
324   PetscFunctionReturn(0);
325 }
326 
327 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
328 {
329   PetscBool      isascii;
330 
331   PetscFunctionBegin;
332   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
333   if (isascii) {
334     TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
335     PetscCall(PetscViewerASCIIPrintf(viewer,"  frtol - relative tolerance in function value %g\n",(double)pseudo->frtol));
336     PetscCall(PetscViewerASCIIPrintf(viewer,"  fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol));
337     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_initial - initial timestep %g\n",(double)pseudo->dt_initial));
338     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment));
339     PetscCall(PetscViewerASCIIPrintf(viewer,"  dt_max - maximum time %g\n",(double)pseudo->dt_max));
340   }
341   PetscFunctionReturn(0);
342 }
343 
344 /* ----------------------------------------------------------------------------- */
345 /*@C
346    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
347    last timestep.
348 
349    Logically Collective on TS
350 
351    Input Parameters:
352 +  ts - timestep context
353 .  dt - user-defined function to verify timestep
354 -  ctx - [optional] user-defined context for private data
355          for the timestep verification routine (may be NULL)
356 
357    Level: advanced
358 
359    Calling sequence of func:
360 $  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
361 
362 +  update - latest solution vector
363 .  ctx - [optional] timestep context
364 .  newdt - the timestep to use for the next step
365 -  flag - flag indicating whether the last time step was acceptable
366 
367    Notes:
368    The routine set here will be called by TSPseudoVerifyTimeStep()
369    during the timestepping process.
370 
371 .seealso: `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()`
372 @*/
373 PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
374 {
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
377   PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));
378   PetscFunctionReturn(0);
379 }
380 
381 /*@
382     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
383     dt when using the TSPseudoTimeStepDefault() routine.
384 
385    Logically Collective on TS
386 
387     Input Parameters:
388 +   ts - the timestep context
389 -   inc - the scaling factor >= 1.0
390 
391     Options Database Key:
392 .    -ts_pseudo_increment <increment> - set pseudo increment
393 
394     Level: advanced
395 
396 .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
397 @*/
398 PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
399 {
400   PetscFunctionBegin;
401   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
402   PetscValidLogicalCollectiveReal(ts,inc,2);
403   PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));
404   PetscFunctionReturn(0);
405 }
406 
407 /*@
408     TSPseudoSetMaxTimeStep - Sets the maximum time step
409     when using the TSPseudoTimeStepDefault() routine.
410 
411    Logically Collective on TS
412 
413     Input Parameters:
414 +   ts - the timestep context
415 -   maxdt - the maximum time step, use a non-positive value to deactivate
416 
417     Options Database Key:
418 .    -ts_pseudo_max_dt <increment> - set pseudo max dt
419 
420     Level: advanced
421 
422 .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
423 @*/
424 PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
425 {
426   PetscFunctionBegin;
427   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
428   PetscValidLogicalCollectiveReal(ts,maxdt,2);
429   PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
430   PetscFunctionReturn(0);
431 }
432 
433 /*@
434     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
435     is computed via the formula
436 $         dt = initial_dt*initial_fnorm/current_fnorm
437       rather than the default update,
438 $         dt = current_dt*previous_fnorm/current_fnorm.
439 
440    Logically Collective on TS
441 
442     Input Parameter:
443 .   ts - the timestep context
444 
445     Options Database Key:
446 .    -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment
447 
448     Level: advanced
449 
450 .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()`
451 @*/
452 PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
453 {
454   PetscFunctionBegin;
455   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
456   PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));
457   PetscFunctionReturn(0);
458 }
459 
460 /*@C
461    TSPseudoSetTimeStep - Sets the user-defined routine to be
462    called at each pseudo-timestep to update the timestep.
463 
464    Logically Collective on TS
465 
466    Input Parameters:
467 +  ts - timestep context
468 .  dt - function to compute timestep
469 -  ctx - [optional] user-defined context for private data
470          required by the function (may be NULL)
471 
472    Level: intermediate
473 
474    Calling sequence of func:
475 $  func (TS ts,PetscReal *newdt,void *ctx);
476 
477 +  newdt - the newly computed timestep
478 -  ctx - [optional] timestep context
479 
480    Notes:
481    The routine set here will be called by TSPseudoComputeTimeStep()
482    during the timestepping process.
483    If not set then TSPseudoTimeStepDefault() is automatically used
484 
485 .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()`
486 @*/
487 PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
488 {
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
491   PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));
492   PetscFunctionReturn(0);
493 }
494 
495 /* ----------------------------------------------------------------------------- */
496 
497 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
498 static PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
499 {
500   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
501 
502   PetscFunctionBegin;
503   pseudo->verify    = dt;
504   pseudo->verifyctx = ctx;
505   PetscFunctionReturn(0);
506 }
507 
508 static PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
509 {
510   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
511 
512   PetscFunctionBegin;
513   pseudo->dt_increment = inc;
514   PetscFunctionReturn(0);
515 }
516 
517 static PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
518 {
519   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
520 
521   PetscFunctionBegin;
522   pseudo->dt_max = maxdt;
523   PetscFunctionReturn(0);
524 }
525 
526 static PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
527 {
528   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
529 
530   PetscFunctionBegin;
531   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
532   PetscFunctionReturn(0);
533 }
534 
535 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
536 static PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
537 {
538   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
539 
540   PetscFunctionBegin;
541   pseudo->dt    = dt;
542   pseudo->dtctx = ctx;
543   PetscFunctionReturn(0);
544 }
545 
546 /* ----------------------------------------------------------------------------- */
547 /*MC
548       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
549 
550   This method solves equations of the form
551 
552 $    F(X,Xdot) = 0
553 
554   for steady state using the iteration
555 
556 $    [G'] S = -F(X,0)
557 $    X += S
558 
559   where
560 
561 $    G(Y) = F(Y,(Y-X)/dt)
562 
563   This is linearly-implicit Euler with the residual always evaluated "at steady
564   state".  See note below.
565 
566   Options database keys:
567 +  -ts_pseudo_increment <real> - ratio of increase dt
568 .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
569 .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
570 -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
571 
572   Level: beginner
573 
574   References:
575 +  * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
576 -  * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
577 
578   Notes:
579   The residual computed by this method includes the transient term (Xdot is computed instead of
580   always being zero), but since the prediction from the last step is always the solution from the
581   last step, on the first Newton iteration we have
582 
583 $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
584 
585   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
586   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
587   algorithm is no longer the one described in the referenced papers.
588 
589 .seealso: `TSCreate()`, `TS`, `TSSetType()`
590 
591 M*/
592 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
593 {
594   TS_Pseudo      *pseudo;
595   SNES           snes;
596   SNESType       stype;
597 
598   PetscFunctionBegin;
599   ts->ops->reset          = TSReset_Pseudo;
600   ts->ops->destroy        = TSDestroy_Pseudo;
601   ts->ops->view           = TSView_Pseudo;
602   ts->ops->setup          = TSSetUp_Pseudo;
603   ts->ops->step           = TSStep_Pseudo;
604   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
605   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
606   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
607   ts->default_adapt_type  = TSADAPTNONE;
608   ts->usessnes            = PETSC_TRUE;
609 
610   PetscCall(TSGetSNES(ts,&snes));
611   PetscCall(SNESGetType(snes,&stype));
612   if (!stype) PetscCall(SNESSetType(snes,SNESKSPONLY));
613 
614   PetscCall(PetscNewLog(ts,&pseudo));
615   ts->data = (void*)pseudo;
616 
617   pseudo->dt                           = TSPseudoTimeStepDefault;
618   pseudo->dtctx                        = NULL;
619   pseudo->dt_increment                 = 1.1;
620   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
621   pseudo->fnorm                        = -1;
622   pseudo->fnorm_initial                = -1;
623   pseudo->fnorm_previous               = -1;
624  #if defined(PETSC_USE_REAL_SINGLE)
625   pseudo->fatol                        = 1.e-25;
626   pseudo->frtol                        = 1.e-5;
627 #else
628   pseudo->fatol                        = 1.e-50;
629   pseudo->frtol                        = 1.e-12;
630 #endif
631   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo));
632   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo));
633   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo));
634   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo));
635   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo));
636   PetscFunctionReturn(0);
637 }
638 
639 /*@C
640    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
641    Use with TSPseudoSetTimeStep().
642 
643    Collective on TS
644 
645    Input Parameters:
646 +  ts - the timestep context
647 -  dtctx - unused timestep context
648 
649    Output Parameter:
650 .  newdt - the timestep to use for the next step
651 
652    Level: advanced
653 
654 .seealso: `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`
655 @*/
656 PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
657 {
658   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
659   PetscReal      inc = pseudo->dt_increment;
660 
661   PetscFunctionBegin;
662   PetscCall(VecZeroEntries(pseudo->xdot));
663   PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE));
664   PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm));
665   if (pseudo->fnorm_initial < 0) {
666     /* first time through so compute initial function norm */
667     pseudo->fnorm_initial  = pseudo->fnorm;
668     pseudo->fnorm_previous = pseudo->fnorm;
669   }
670   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
671   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
672   else                                           *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm;
673   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
674   pseudo->fnorm_previous = pseudo->fnorm;
675   PetscFunctionReturn(0);
676 }
677