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