xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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   PetscReal fatol,frtol;
26 } TS_Pseudo;
27 
28 /* ------------------------------------------------------------------------------*/
29 
30 /*@C
31     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32     pseudo-timestepping process.
33 
34     Collective on TS
35 
36     Input Parameter:
37 .   ts - timestep context
38 
39     Output Parameter:
40 .   dt - newly computed timestep
41 
42     Level: developer
43 
44     Notes:
45     The routine to be called here to compute the timestep should be
46     set by calling TSPseudoSetTimeStep().
47 
48 .keywords: timestep, pseudo, compute
49 
50 .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep()
51 @*/
52 PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
53 {
54   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
59   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
60   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
61   PetscFunctionReturn(0);
62 }
63 
64 
65 /* ------------------------------------------------------------------------------*/
66 /*@C
67    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
68 
69    Collective on TS
70 
71    Input Parameters:
72 +  ts - the timestep context
73 .  dtctx - unused timestep context
74 -  update - latest solution vector
75 
76    Output Parameters:
77 +  newdt - the timestep to use for the next step
78 -  flag - flag indicating whether the last time step was acceptable
79 
80    Level: advanced
81 
82    Note:
83    This routine always returns a flag of 1, indicating an acceptable
84    timestep.
85 
86 .keywords: timestep, pseudo, default, verify
87 
88 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
89 @*/
90 PetscErrorCode  TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
91 {
92   PetscFunctionBegin;
93   *flag = PETSC_TRUE;
94   PetscFunctionReturn(0);
95 }
96 
97 
98 /*@
99     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
100 
101     Collective on TS
102 
103     Input Parameters:
104 +   ts - timestep context
105 -   update - latest solution vector
106 
107     Output Parameters:
108 +   dt - newly computed timestep (if it had to shrink)
109 -   flag - indicates if current timestep was ok
110 
111     Level: advanced
112 
113     Notes:
114     The routine to be called here to compute the timestep should be
115     set by calling TSPseudoSetVerifyTimeStep().
116 
117 .keywords: timestep, pseudo, verify
118 
119 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
120 @*/
121 PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag)
122 {
123   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
124   PetscErrorCode ierr;
125 
126   PetscFunctionBegin;
127   *flag = PETSC_TRUE;
128   if(pseudo->verify) {
129     ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
130   }
131   PetscFunctionReturn(0);
132 }
133 
134 /* --------------------------------------------------------------------------------*/
135 
136 static PetscErrorCode TSStep_Pseudo(TS ts)
137 {
138   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
139   PetscInt            nits,lits,reject;
140   PetscBool           stepok;
141   PetscReal           next_time_step = ts->time_step;
142   PetscErrorCode      ierr;
143 
144   PetscFunctionBegin;
145   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
146   ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr);
147   ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr);
148   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
149     ts->time_step = next_time_step;
150     ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr);
151     ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr);
152     ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
153     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
154     ts->snes_its += nits; ts->ksp_its += lits;
155     ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr);
156     ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok);CHKERRQ(ierr);
157     if (!stepok) {next_time_step = ts->time_step; continue;}
158     pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
159     ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr);
160     if (stepok) break;
161   }
162   if (reject >= ts->max_reject) {
163     ts->reason = TS_DIVERGED_STEP_REJECTED;
164     ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
165     PetscFunctionReturn(0);
166   }
167 
168   ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr);
169   ts->ptime += ts->time_step;
170   ts->time_step = next_time_step;
171 
172   if (pseudo->fnorm < 0) {
173     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
174     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
175     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
176   }
177   if (pseudo->fnorm < pseudo->fatol) {
178     ts->reason = TS_CONVERGED_PSEUDO_FATOL;
179     ierr = PetscInfo3(ts,"Step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr);
180     PetscFunctionReturn(0);
181   }
182   if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) {
183     ts->reason = TS_CONVERGED_PSEUDO_FRTOL;
184     ierr = PetscInfo4(ts,"Step=%D, converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,pseudo->fnorm,pseudo->fnorm_initial,pseudo->fatol);CHKERRQ(ierr);
185     PetscFunctionReturn(0);
186   }
187   PetscFunctionReturn(0);
188 }
189 
190 /*------------------------------------------------------------*/
191 static PetscErrorCode TSReset_Pseudo(TS ts)
192 {
193   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
194   PetscErrorCode ierr;
195 
196   PetscFunctionBegin;
197   ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr);
198   ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr);
199   ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
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   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr);
211   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr);
212   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
213   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr);
214   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 /*------------------------------------------------------------*/
219 
220 /*
221     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
222 */
223 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
224 {
225   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
226   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
227   PetscScalar       *xdot;
228   PetscErrorCode    ierr;
229   PetscInt          i,n;
230 
231   PetscFunctionBegin;
232   ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
233   ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr);
234   ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
235   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
236   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
237   ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr);
238   ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr);
239   ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr);
240   *Xdot = pseudo->xdot;
241   PetscFunctionReturn(0);
242 }
243 
244 /*
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+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 
272 /*
273    This constructs the Jacobian needed for SNES.  For DAE, this is
274 
275        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
276 
277     and for ODE:
278 
279        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
280 */
281 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
282 {
283   Vec            Xdot;
284   PetscErrorCode ierr;
285 
286   PetscFunctionBegin;
287   ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr);
288   ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr);
289   PetscFunctionReturn(0);
290 }
291 
292 
293 static PetscErrorCode TSSetUp_Pseudo(TS ts)
294 {
295   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
300   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
301   ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 /*------------------------------------------------------------*/
305 
306 static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
307 {
308   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
309   PetscErrorCode ierr;
310   PetscViewer    viewer = (PetscViewer) dummy;
311 
312   PetscFunctionBegin;
313   if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
314     ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
315     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
316     ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
317   }
318   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
319   ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);CHKERRQ(ierr);
320   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323 
324 static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts)
325 {
326   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
327   PetscErrorCode ierr;
328   PetscBool      flg = PETSC_FALSE;
329   PetscViewer    viewer;
330 
331   PetscFunctionBegin;
332   ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr);
333   ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL);CHKERRQ(ierr);
334   if (flg) {
335     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr);
336     ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
337   }
338   flg  = pseudo->increment_dt_from_initial_dt;
339   ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr);
340   pseudo->increment_dt_from_initial_dt = flg;
341   ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr);
342   ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr);
343   ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr);
344   ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr);
345   ierr = PetscOptionsTail();CHKERRQ(ierr);
346   PetscFunctionReturn(0);
347 }
348 
349 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
350 {
351   PetscErrorCode ierr;
352   PetscBool      isascii;
353 
354   PetscFunctionBegin;
355   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
356   if (isascii) {
357     TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
358     ierr = PetscViewerASCIIPrintf(viewer,"Parameters for pseudo timestepping\n");CHKERRQ(ierr);
359     ierr = PetscViewerASCIIPrintf(viewer,"  frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr);
360     ierr = PetscViewerASCIIPrintf(viewer,"  fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr);
361     ierr = PetscViewerASCIIPrintf(viewer,"  dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr);
362     ierr = PetscViewerASCIIPrintf(viewer,"  dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr);
363     ierr = PetscViewerASCIIPrintf(viewer,"  dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr);
364   }
365   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
366   PetscFunctionReturn(0);
367 }
368 
369 /* ----------------------------------------------------------------------------- */
370 /*@C
371    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
372    last timestep.
373 
374    Logically Collective on TS
375 
376    Input Parameters:
377 +  ts - timestep context
378 .  dt - user-defined function to verify timestep
379 -  ctx - [optional] user-defined context for private data
380          for the timestep verification routine (may be NULL)
381 
382    Level: advanced
383 
384    Calling sequence of func:
385 .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);
386 
387 .  update - latest solution vector
388 .  ctx - [optional] timestep context
389 .  newdt - the timestep to use for the next step
390 .  flag - flag indicating whether the last time step was acceptable
391 
392    Notes:
393    The routine set here will be called by TSPseudoVerifyTimeStep()
394    during the timestepping process.
395 
396 .keywords: timestep, pseudo, set, verify
397 
398 .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
399 @*/
400 PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
401 {
402   PetscErrorCode ierr;
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
406   ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr);
407   PetscFunctionReturn(0);
408 }
409 
410 /*@
411     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
412     dt when using the TSPseudoTimeStepDefault() routine.
413 
414    Logically Collective on TS
415 
416     Input Parameters:
417 +   ts - the timestep context
418 -   inc - the scaling factor >= 1.0
419 
420     Options Database Key:
421 .    -ts_pseudo_increment <increment>
422 
423     Level: advanced
424 
425 .keywords: timestep, pseudo, set, increment
426 
427 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
428 @*/
429 PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
430 {
431   PetscErrorCode ierr;
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
435   PetscValidLogicalCollectiveReal(ts,inc,2);
436   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr);
437   PetscFunctionReturn(0);
438 }
439 
440 /*@
441     TSPseudoSetMaxTimeStep - Sets the maximum time step
442     when using the TSPseudoTimeStepDefault() routine.
443 
444    Logically Collective on TS
445 
446     Input Parameters:
447 +   ts - the timestep context
448 -   maxdt - the maximum time step, use a non-positive value to deactivate
449 
450     Options Database Key:
451 .    -ts_pseudo_max_dt <increment>
452 
453     Level: advanced
454 
455 .keywords: timestep, pseudo, set
456 
457 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
458 @*/
459 PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
460 {
461   PetscErrorCode ierr;
462 
463   PetscFunctionBegin;
464   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
465   PetscValidLogicalCollectiveReal(ts,maxdt,2);
466   ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
467   PetscFunctionReturn(0);
468 }
469 
470 /*@
471     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
472     is computed via the formula
473 $         dt = initial_dt*initial_fnorm/current_fnorm
474       rather than the default update,
475 $         dt = current_dt*previous_fnorm/current_fnorm.
476 
477    Logically Collective on TS
478 
479     Input Parameter:
480 .   ts - the timestep context
481 
482     Options Database Key:
483 .    -ts_pseudo_increment_dt_from_initial_dt
484 
485     Level: advanced
486 
487 .keywords: timestep, pseudo, set, increment
488 
489 .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
490 @*/
491 PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
492 {
493   PetscErrorCode ierr;
494 
495   PetscFunctionBegin;
496   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
497   ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 
502 /*@C
503    TSPseudoSetTimeStep - Sets the user-defined routine to be
504    called at each pseudo-timestep to update the timestep.
505 
506    Logically Collective on TS
507 
508    Input Parameters:
509 +  ts - timestep context
510 .  dt - function to compute timestep
511 -  ctx - [optional] user-defined context for private data
512          required by the function (may be NULL)
513 
514    Level: intermediate
515 
516    Calling sequence of func:
517 .  func (TS ts,PetscReal *newdt,void *ctx);
518 
519 .  newdt - the newly computed timestep
520 .  ctx - [optional] timestep context
521 
522    Notes:
523    The routine set here will be called by TSPseudoComputeTimeStep()
524    during the timestepping process.
525    If not set then TSPseudoTimeStepDefault() is automatically used
526 
527 .keywords: timestep, pseudo, set
528 
529 .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
530 @*/
531 PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
532 {
533   PetscErrorCode ierr;
534 
535   PetscFunctionBegin;
536   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
537   ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr);
538   PetscFunctionReturn(0);
539 }
540 
541 /* ----------------------------------------------------------------------------- */
542 
543 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
544 static PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
545 {
546   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
547 
548   PetscFunctionBegin;
549   pseudo->verify    = dt;
550   pseudo->verifyctx = ctx;
551   PetscFunctionReturn(0);
552 }
553 
554 static PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
555 {
556   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
557 
558   PetscFunctionBegin;
559   pseudo->dt_increment = inc;
560   PetscFunctionReturn(0);
561 }
562 
563 static PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
564 {
565   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
566 
567   PetscFunctionBegin;
568   pseudo->dt_max = maxdt;
569   PetscFunctionReturn(0);
570 }
571 
572 static PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
573 {
574   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
575 
576   PetscFunctionBegin;
577   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
578   PetscFunctionReturn(0);
579 }
580 
581 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
582 static PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
583 {
584   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
585 
586   PetscFunctionBegin;
587   pseudo->dt    = dt;
588   pseudo->dtctx = ctx;
589   PetscFunctionReturn(0);
590 }
591 
592 /* ----------------------------------------------------------------------------- */
593 /*MC
594       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
595 
596   This method solves equations of the form
597 
598 $    F(X,Xdot) = 0
599 
600   for steady state using the iteration
601 
602 $    [G'] S = -F(X,0)
603 $    X += S
604 
605   where
606 
607 $    G(Y) = F(Y,(Y-X)/dt)
608 
609   This is linearly-implicit Euler with the residual always evaluated "at steady
610   state".  See note below.
611 
612   Options database keys:
613 +  -ts_pseudo_increment <real> - ratio of increase dt
614 .  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
615 .  -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
616 -  -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol
617 
618   Level: beginner
619 
620   References:
621 +  1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
622 -  2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
623 
624   Notes:
625   The residual computed by this method includes the transient term (Xdot is computed instead of
626   always being zero), but since the prediction from the last step is always the solution from the
627   last step, on the first Newton iteration we have
628 
629 $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
630 
631   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
632   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
633   algorithm is no longer the one described in the referenced papers.
634 
635 .seealso:  TSCreate(), TS, TSSetType()
636 
637 M*/
638 PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
639 {
640   TS_Pseudo      *pseudo;
641   PetscErrorCode ierr;
642   SNES           snes;
643   SNESType       stype;
644 
645   PetscFunctionBegin;
646   ts->ops->reset          = TSReset_Pseudo;
647   ts->ops->destroy        = TSDestroy_Pseudo;
648   ts->ops->view           = TSView_Pseudo;
649   ts->ops->setup          = TSSetUp_Pseudo;
650   ts->ops->step           = TSStep_Pseudo;
651   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
652   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
653   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;
654 
655   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
656   ierr = SNESGetType(snes,&stype);CHKERRQ(ierr);
657   if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
658 
659   ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr);
660   ts->data = (void*)pseudo;
661 
662   pseudo->dt                           = TSPseudoTimeStepDefault;
663   pseudo->dtctx                        = NULL;
664   pseudo->dt_increment                 = 1.1;
665   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
666   pseudo->fnorm                        = -1;
667   pseudo->fnorm_initial                = -1;
668   pseudo->fnorm_previous               = -1;
669  #if defined(PETSC_USE_REAL_SINGLE)
670   pseudo->fatol                        = 1.e-25;
671   pseudo->frtol                        = 1.e-5;
672 #else
673   pseudo->fatol                        = 1.e-50;
674   pseudo->frtol                        = 1.e-12;
675 #endif
676   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
677   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
678   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr);
679   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
680   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
681   PetscFunctionReturn(0);
682 }
683 
684 /*@C
685    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
686    Use with TSPseudoSetTimeStep().
687 
688    Collective on TS
689 
690    Input Parameters:
691 .  ts - the timestep context
692 .  dtctx - unused timestep context
693 
694    Output Parameter:
695 .  newdt - the timestep to use for the next step
696 
697    Level: advanced
698 
699 .keywords: timestep, pseudo, default
700 
701 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
702 @*/
703 PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
704 {
705   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
706   PetscReal      inc = pseudo->dt_increment;
707   PetscErrorCode ierr;
708 
709   PetscFunctionBegin;
710   ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr);
711   ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr);
712   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
713   if (pseudo->fnorm_initial < 0) {
714     /* first time through so compute initial function norm */
715     pseudo->fnorm_initial  = pseudo->fnorm;
716     pseudo->fnorm_previous = pseudo->fnorm;
717   }
718   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
719   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
720   else                                           *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm;
721   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
722   pseudo->fnorm_previous = pseudo->fnorm;
723   PetscFunctionReturn(0);
724 }
725