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