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