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