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