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