xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 7ba3a57c1a41d0d4ac9568ddc4edec2fa54aa2cd)
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   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
208   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr);
209   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr);
210   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 /*------------------------------------------------------------*/
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "TSPseudoGetXdot"
218 /*
219     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
220 */
221 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
222 {
223   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
224   const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
225   PetscScalar    *xdot;
226   PetscErrorCode ierr;
227   PetscInt       i,n;
228 
229   PetscFunctionBegin;
230   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
231   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
232   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
233   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
234   for (i=0; i<n; i++) {
235     xdot[i] = mdt*(xnp1[i] - xn[i]);
236   }
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 #undef __FUNCT__
245 #define __FUNCT__ "SNESTSFormFunction_Pseudo"
246 /*
247     The transient residual is
248 
249         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
250 
251     or for ODE,
252 
253         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
254 
255     This is the function that must be evaluated for transient simulation and for
256     finite difference Jacobians.  On the first Newton step, this algorithm uses
257     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
258     residual is actually the steady state residual.  Pseudotransient
259     continuation as described in the literature is a linearly implicit
260     algorithm, it only takes this one Newton step with the steady state
261     residual, and then advances to the next time step.
262 */
263 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
264 {
265   Vec            Xdot;
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
270   ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
271   PetscFunctionReturn(0);
272 }
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "SNESTSFormJacobian_Pseudo"
276 /*
277    This constructs the Jacobian needed for SNES.  For DAE, this is
278 
279        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
280 
281     and for ODE:
282 
283        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
284 */
285 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts)
286 {
287   Vec            Xdot;
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
292   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "TSSetUp_Pseudo"
299 static PetscErrorCode TSSetUp_Pseudo(TS ts)
300 {
301   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
306   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
307   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
308   PetscFunctionReturn(0);
309 }
310 /*------------------------------------------------------------*/
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "TSPseudoMonitorDefault"
314 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
315 {
316   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
317   PetscErrorCode   ierr;
318   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
319 
320   PetscFunctionBegin;
321   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
322     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
323     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
324     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
325   }
326   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
327   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
328   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "TSSetFromOptions_Pseudo"
334 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
335 {
336   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
337   PetscErrorCode  ierr;
338   PetscBool       flg = PETSC_FALSE;
339   PetscViewer     viewer;
340 
341   PetscFunctionBegin;
342   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
343     ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
344     if (flg) {
345       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr);
346       ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
347     }
348     flg  = PETSC_FALSE;
349     ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
350     if (flg) {
351       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
352     }
353     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
354 
355     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
356   ierr = PetscOptionsTail();CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "TSView_Pseudo"
362 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
363 {
364   PetscErrorCode ierr;
365 
366   PetscFunctionBegin;
367   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
368   PetscFunctionReturn(0);
369 }
370 
371 /* ----------------------------------------------------------------------------- */
372 #undef __FUNCT__
373 #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
374 /*@C
375    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
376    last timestep.
377 
378    Logically Collective on TS
379 
380    Input Parameters:
381 +  ts - timestep context
382 .  dt - user-defined function to verify timestep
383 -  ctx - [optional] user-defined context for private data
384          for the timestep verification routine (may be PETSC_NULL)
385 
386    Level: advanced
387 
388    Calling sequence of func:
389 .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
390 
391 .  update - latest solution vector
392 .  ctx - [optional] timestep context
393 .  newdt - the timestep to use for the next step
394 .  flag - flag indicating whether the last time step was acceptable
395 
396    Notes:
397    The routine set here will be called by TSPseudoVerifyTimeStep()
398    during the timestepping process.
399 
400 .keywords: timestep, pseudo, set, verify
401 
402 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
403 @*/
404 PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx)
405 {
406   PetscErrorCode ierr;
407 
408   PetscFunctionBegin;
409   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
410   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool  *),void *),(ts,dt,ctx));CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
416 /*@
417     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
418     dt when using the TSPseudoDefaultTimeStep() routine.
419 
420    Logically Collective on TS
421 
422     Input Parameters:
423 +   ts - the timestep context
424 -   inc - the scaling factor >= 1.0
425 
426     Options Database Key:
427 $    -ts_pseudo_increment <increment>
428 
429     Level: advanced
430 
431 .keywords: timestep, pseudo, set, increment
432 
433 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
434 @*/
435 PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
436 {
437   PetscErrorCode ierr;
438 
439   PetscFunctionBegin;
440   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
441   PetscValidLogicalCollectiveReal(ts,inc,2);
442   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
443   PetscFunctionReturn(0);
444 }
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
448 /*@
449     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
450     is computed via the formula
451 $         dt = initial_dt*initial_fnorm/current_fnorm
452       rather than the default update,
453 $         dt = current_dt*previous_fnorm/current_fnorm.
454 
455    Logically Collective on TS
456 
457     Input Parameter:
458 .   ts - the timestep context
459 
460     Options Database Key:
461 $    -ts_pseudo_increment_dt_from_initial_dt
462 
463     Level: advanced
464 
465 .keywords: timestep, pseudo, set, increment
466 
467 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
468 @*/
469 PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
470 {
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
475   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
476   PetscFunctionReturn(0);
477 }
478 
479 
480 #undef __FUNCT__
481 #define __FUNCT__ "TSPseudoSetTimeStep"
482 /*@C
483    TSPseudoSetTimeStep - Sets the user-defined routine to be
484    called at each pseudo-timestep to update the timestep.
485 
486    Logically Collective on TS
487 
488    Input Parameters:
489 +  ts - timestep context
490 .  dt - function to compute timestep
491 -  ctx - [optional] user-defined context for private data
492          required by the function (may be PETSC_NULL)
493 
494    Level: intermediate
495 
496    Calling sequence of func:
497 .  func (TS ts,PetscReal *newdt,void *ctx);
498 
499 .  newdt - the newly computed timestep
500 .  ctx - [optional] timestep context
501 
502    Notes:
503    The routine set here will be called by TSPseudoComputeTimeStep()
504    during the timestepping process.
505 
506 .keywords: timestep, pseudo, set
507 
508 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
509 @*/
510 PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
511 {
512   PetscErrorCode ierr;
513 
514   PetscFunctionBegin;
515   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
516   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr);
517   PetscFunctionReturn(0);
518 }
519 
520 /* ----------------------------------------------------------------------------- */
521 
522 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
523 EXTERN_C_BEGIN
524 #undef __FUNCT__
525 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
526 PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
527 {
528   TS_Pseudo *pseudo;
529 
530   PetscFunctionBegin;
531   pseudo              = (TS_Pseudo*)ts->data;
532   pseudo->verify      = dt;
533   pseudo->verifyctx   = ctx;
534   PetscFunctionReturn(0);
535 }
536 EXTERN_C_END
537 
538 EXTERN_C_BEGIN
539 #undef __FUNCT__
540 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
541 PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
542 {
543   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
544 
545   PetscFunctionBegin;
546   pseudo->dt_increment = inc;
547   PetscFunctionReturn(0);
548 }
549 EXTERN_C_END
550 
551 EXTERN_C_BEGIN
552 #undef __FUNCT__
553 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
554 PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
555 {
556   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
557 
558   PetscFunctionBegin;
559   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
560   PetscFunctionReturn(0);
561 }
562 EXTERN_C_END
563 
564 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
565 EXTERN_C_BEGIN
566 #undef __FUNCT__
567 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
568 PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
569 {
570   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
571 
572   PetscFunctionBegin;
573   pseudo->dt      = dt;
574   pseudo->dtctx   = ctx;
575   PetscFunctionReturn(0);
576 }
577 EXTERN_C_END
578 
579 /* ----------------------------------------------------------------------------- */
580 /*MC
581       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
582 
583   This method solves equations of the form
584 
585 $    F(X,Xdot) = 0
586 
587   for steady state using the iteration
588 
589 $    [G'] S = -F(X,0)
590 $    X += S
591 
592   where
593 
594 $    G(Y) = F(Y,(Y-X)/dt)
595 
596   This is linearly-implicit Euler with the residual always evaluated "at steady
597   state".  See note below.
598 
599   Options database keys:
600 +  -ts_pseudo_increment <real> - ratio of increase dt
601 -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
602 
603   Level: beginner
604 
605   References:
606   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
607   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
608 
609   Notes:
610   The residual computed by this method includes the transient term (Xdot is computed instead of
611   always being zero), but since the prediction from the last step is always the solution from the
612   last step, on the first Newton iteration we have
613 
614 $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
615 
616   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
617   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
618   algorithm is no longer the one described in the referenced papers.
619 
620 .seealso:  TSCreate(), TS, TSSetType()
621 
622 M*/
623 EXTERN_C_BEGIN
624 #undef __FUNCT__
625 #define __FUNCT__ "TSCreate_Pseudo"
626 PetscErrorCode  TSCreate_Pseudo(TS ts)
627 {
628   TS_Pseudo      *pseudo;
629   PetscErrorCode ierr;
630   SNES           snes;
631   const SNESType stype;
632 
633   PetscFunctionBegin;
634   ts->ops->reset           = TSReset_Pseudo;
635   ts->ops->destroy         = TSDestroy_Pseudo;
636   ts->ops->view            = TSView_Pseudo;
637 
638   ts->ops->setup           = TSSetUp_Pseudo;
639   ts->ops->step            = TSStep_Pseudo;
640   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
641   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
642   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;
643 
644   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
645   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
646   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
647 
648   ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr);
649   ts->data = (void*)pseudo;
650 
651   pseudo->dt_increment                 = 1.1;
652   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
653   pseudo->dt                           = TSPseudoDefaultTimeStep;
654   pseudo->fnorm                        = -1;
655 
656   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
657                     "TSPseudoSetVerifyTimeStep_Pseudo",
658                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
659   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
660                     "TSPseudoSetTimeStepIncrement_Pseudo",
661                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
662   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
663                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
664                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
665   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
666                     "TSPseudoSetTimeStep_Pseudo",
667                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
668   PetscFunctionReturn(0);
669 }
670 EXTERN_C_END
671 
672 #undef __FUNCT__
673 #define __FUNCT__ "TSPseudoDefaultTimeStep"
674 /*@C
675    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
676    Use with TSPseudoSetTimeStep().
677 
678    Collective on TS
679 
680    Input Parameters:
681 .  ts - the timestep context
682 .  dtctx - unused timestep context
683 
684    Output Parameter:
685 .  newdt - the timestep to use for the next step
686 
687    Level: advanced
688 
689 .keywords: timestep, pseudo, default
690 
691 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
692 @*/
693 PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
694 {
695   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
696   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
697   PetscErrorCode ierr;
698 
699   PetscFunctionBegin;
700   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
701   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
702   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
703   if (pseudo->initial_fnorm == 0.0) {
704     /* first time through so compute initial function norm */
705     pseudo->initial_fnorm = pseudo->fnorm;
706     fnorm_previous        = pseudo->fnorm;
707   }
708   if (pseudo->fnorm == 0.0) {
709     *newdt = 1.e12*inc*ts->time_step;
710   } else if (pseudo->increment_dt_from_initial_dt) {
711     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
712   } else {
713     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
714   }
715   pseudo->fnorm_previous = pseudo->fnorm;
716   PetscFunctionReturn(0);
717 }
718