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