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