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