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