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