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