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