xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 4482741e5b2e2bbc854fb1f8dba65221386520f2)
1 /*$Id: posindep.c,v 1.60 2001/09/11 16:34:22 bsmith Exp $*/
2 /*
3        Code for Timestepping with implicit backwards Euler.
4 */
5 #include "src/ts/tsimpl.h"                /*I   "petscts.h"   I*/
6 
7 typedef struct {
8   Vec  update;      /* work vector where new solution is formed */
9   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
10   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
11 
12   /* information used for Pseudo-timestepping */
13 
14   int    (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
15   void   *dtctx;
16   int    (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */
17   void   *verifyctx;
18 
19   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
20   PetscReal  fnorm_previous;
21 
22   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
23   PetscTruth increment_dt_from_initial_dt;
24 } TS_Pseudo;
25 
26 /* ------------------------------------------------------------------------------*/
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "TSPseudoComputeTimeStep"
30 /*@
31     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
32     pseudo-timestepping process.
33 
34     Collective on TS
35 
36     Input Parameter:
37 .   ts - timestep context
38 
39     Output Parameter:
40 .   dt - newly computed timestep
41 
42     Level: advanced
43 
44     Notes:
45     The routine to be called here to compute the timestep should be
46     set by calling TSPseudoSetTimeStep().
47 
48 .keywords: timestep, pseudo, compute
49 
50 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
51 @*/
52 int TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
53 {
54   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
55   int       ierr;
56 
57   PetscFunctionBegin;
58   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
59   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
60   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
61   PetscFunctionReturn(0);
62 }
63 
64 
65 /* ------------------------------------------------------------------------------*/
66 #undef __FUNCT__
67 #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
68 /*@C
69    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
70 
71    Collective on TS
72 
73    Input Parameters:
74 +  ts - the timestep context
75 .  dtctx - unused timestep context
76 -  update - latest solution vector
77 
78    Output Parameters:
79 +  newdt - the timestep to use for the next step
80 -  flag - flag indicating whether the last time step was acceptable
81 
82    Level: advanced
83 
84    Note:
85    This routine always returns a flag of 1, indicating an acceptable
86    timestep.
87 
88 .keywords: timestep, pseudo, default, verify
89 
90 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
91 @*/
92 int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag)
93 {
94   PetscFunctionBegin;
95   *flag = 1;
96   PetscFunctionReturn(0);
97 }
98 
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "TSPseudoVerifyTimeStep"
102 /*@
103     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
104 
105     Collective on TS
106 
107     Input Parameters:
108 +   ts - timestep context
109 -   update - latest solution vector
110 
111     Output Parameters:
112 +   dt - newly computed timestep (if it had to shrink)
113 -   flag - indicates if current timestep was ok
114 
115     Level: advanced
116 
117     Notes:
118     The routine to be called here to compute the timestep should be
119     set by calling TSPseudoSetVerifyTimeStep().
120 
121 .keywords: timestep, pseudo, verify
122 
123 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
124 @*/
125 int TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag)
126 {
127   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
128   int       ierr;
129 
130   PetscFunctionBegin;
131   if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);}
132 
133   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
134 
135   PetscFunctionReturn(0);
136 }
137 
138 /* --------------------------------------------------------------------------------*/
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "TSStep_Pseudo"
142 static int TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime)
143 {
144   Vec       sol = ts->vec_sol;
145   int       ierr,i,max_steps = ts->max_steps,its,ok,lits;
146   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
147   PetscReal current_time_step;
148 
149   PetscFunctionBegin;
150   *steps = -ts->steps;
151 
152   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
153   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
154     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
155     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
156     current_time_step = ts->time_step;
157     while (PETSC_TRUE) {
158       ts->ptime  += current_time_step;
159       ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr);
160       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
161       ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
162       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
163       if (ok) break;
164       ts->ptime        -= current_time_step;
165       current_time_step = ts->time_step;
166     }
167     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
168     ts->steps++;
169   }
170   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
171   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
172   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
173 
174   *steps += ts->steps;
175   *ptime  = ts->ptime;
176   PetscFunctionReturn(0);
177 }
178 
179 /*------------------------------------------------------------*/
180 #undef __FUNCT__
181 #define __FUNCT__ "TSDestroy_Pseudo"
182 static int TSDestroy_Pseudo(TS ts)
183 {
184   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
185   int       ierr;
186 
187   PetscFunctionBegin;
188   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
189   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
190   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
191   ierr = PetscFree(pseudo);CHKERRQ(ierr);
192   PetscFunctionReturn(0);
193 }
194 
195 
196 /*------------------------------------------------------------*/
197 
198 /*
199     This defines the nonlinear equation that is to be solved with SNES
200 
201               (U^{n+1} - U^{n})/dt - F(U^{n+1})
202 */
203 #undef __FUNCT__
204 #define __FUNCT__ "TSPseudoFunction"
205 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
206 {
207   TS     ts = (TS) ctx;
208   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
209   int    ierr,i,n;
210 
211   PetscFunctionBegin;
212   /* apply user provided function */
213   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
214   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
215   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
216   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
217   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
218   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
219   for (i=0; i<n; i++) {
220     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
221   }
222   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
223   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
224   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
225 
226   PetscFunctionReturn(0);
227 }
228 
229 /*
230    This constructs the Jacobian needed for SNES
231 
232              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
233 */
234 #undef __FUNCT__
235 #define __FUNCT__ "TSPseudoJacobian"
236 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
237 {
238   TS            ts = (TS) ctx;
239   int           ierr;
240   PetscScalar   mone = -1.0,mdt = 1.0/ts->time_step;
241 
242   PetscFunctionBegin;
243   /* construct users Jacobian */
244   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
245 
246   /* shift and scale Jacobian */
247   ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
248   ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
249   if (*BB != *AA && *str != SAME_PRECONDITIONER) {
250     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
251     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
252   }
253 
254   PetscFunctionReturn(0);
255 }
256 
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "TSSetUp_Pseudo"
260 static int TSSetUp_Pseudo(TS ts)
261 {
262   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
263   int       ierr;
264 
265   PetscFunctionBegin;
266   /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */
267   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
268   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
269   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
270   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
271   PetscFunctionReturn(0);
272 }
273 /*------------------------------------------------------------*/
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "TSPseudoDefaultMonitor"
277 int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
278 {
279   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
280   int       ierr;
281 
282   PetscFunctionBegin;
283   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "TSSetFromOptions_Pseudo"
289 static int TSSetFromOptions_Pseudo(TS ts)
290 {
291   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
292   int        ierr;
293   PetscTruth flg;
294 
295   PetscFunctionBegin;
296 
297   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
298     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
299     if (flg) {
300       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
301     }
302     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
303     if (flg) {
304       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
305     }
306     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
307   ierr = PetscOptionsTail();CHKERRQ(ierr);
308 
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "TSView_Pseudo"
314 static int TSView_Pseudo(TS ts,PetscViewer viewer)
315 {
316   PetscFunctionBegin;
317   PetscFunctionReturn(0);
318 }
319 
320 /* ----------------------------------------------------------------------------- */
321 #undef __FUNCT__
322 #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
323 /*@C
324    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
325    last timestep.
326 
327    Collective on TS
328 
329    Input Parameters:
330 +  ts - timestep context
331 .  dt - user-defined function to verify timestep
332 -  ctx - [optional] user-defined context for private data
333          for the timestep verification routine (may be PETSC_NULL)
334 
335    Level: advanced
336 
337    Calling sequence of func:
338 .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
339 
340 .  update - latest solution vector
341 .  ctx - [optional] timestep context
342 .  newdt - the timestep to use for the next step
343 .  flag - flag indicating whether the last time step was acceptable
344 
345    Notes:
346    The routine set here will be called by TSPseudoVerifyTimeStep()
347    during the timestepping process.
348 
349 .keywords: timestep, pseudo, set, verify
350 
351 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
352 @*/
353 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
354 {
355   int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *);
356 
357   PetscFunctionBegin;
358   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
359 
360   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
361   if (f) {
362     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
363   }
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
369 /*@
370     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
371     dt when using the TSPseudoDefaultTimeStep() routine.
372 
373    Collective on TS
374 
375     Input Parameters:
376 +   ts - the timestep context
377 -   inc - the scaling factor >= 1.0
378 
379     Options Database Key:
380 $    -ts_pseudo_increment <increment>
381 
382     Level: advanced
383 
384 .keywords: timestep, pseudo, set, increment
385 
386 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
387 @*/
388 int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
389 {
390   int ierr,(*f)(TS,PetscReal);
391 
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
394 
395   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
396   if (f) {
397     ierr = (*f)(ts,inc);CHKERRQ(ierr);
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 #undef __FUNCT__
403 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
404 /*@
405     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
406     is computed via the formula
407 $         dt = initial_dt*initial_fnorm/current_fnorm
408       rather than the default update,
409 $         dt = current_dt*previous_fnorm/current_fnorm.
410 
411    Collective on TS
412 
413     Input Parameter:
414 .   ts - the timestep context
415 
416     Options Database Key:
417 $    -ts_pseudo_increment_dt_from_initial_dt
418 
419     Level: advanced
420 
421 .keywords: timestep, pseudo, set, increment
422 
423 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
424 @*/
425 int TSPseudoIncrementDtFromInitialDt(TS ts)
426 {
427   int ierr,(*f)(TS);
428 
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
431 
432   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
433   if (f) {
434     ierr = (*f)(ts);CHKERRQ(ierr);
435   }
436   PetscFunctionReturn(0);
437 }
438 
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "TSPseudoSetTimeStep"
442 /*@C
443    TSPseudoSetTimeStep - Sets the user-defined routine to be
444    called at each pseudo-timestep to update the timestep.
445 
446    Collective on TS
447 
448    Input Parameters:
449 +  ts - timestep context
450 .  dt - function to compute timestep
451 -  ctx - [optional] user-defined context for private data
452          required by the function (may be PETSC_NULL)
453 
454    Level: intermediate
455 
456    Calling sequence of func:
457 .  func (TS ts,PetscReal *newdt,void *ctx);
458 
459 .  newdt - the newly computed timestep
460 .  ctx - [optional] timestep context
461 
462    Notes:
463    The routine set here will be called by TSPseudoComputeTimeStep()
464    during the timestepping process.
465 
466 .keywords: timestep, pseudo, set
467 
468 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
469 @*/
470 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
471 {
472   int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *);
473 
474   PetscFunctionBegin;
475   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
476 
477   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
478   if (f) {
479     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
480   }
481   PetscFunctionReturn(0);
482 }
483 
484 /* ----------------------------------------------------------------------------- */
485 
486 typedef int (*FCN1)(TS,Vec,void*,PetscReal*,int*); /* force argument to next function to not be extern C*/
487 EXTERN_C_BEGIN
488 #undef __FUNCT__
489 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
490 int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
491 {
492   TS_Pseudo *pseudo;
493 
494   PetscFunctionBegin;
495   pseudo              = (TS_Pseudo*)ts->data;
496   pseudo->verify      = dt;
497   pseudo->verifyctx   = ctx;
498   PetscFunctionReturn(0);
499 }
500 EXTERN_C_END
501 
502 EXTERN_C_BEGIN
503 #undef __FUNCT__
504 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
505 int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
506 {
507   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
508 
509   PetscFunctionBegin;
510   pseudo->dt_increment = inc;
511   PetscFunctionReturn(0);
512 }
513 EXTERN_C_END
514 
515 EXTERN_C_BEGIN
516 #undef __FUNCT__
517 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
518 int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
519 {
520   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
521 
522   PetscFunctionBegin;
523   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
524   PetscFunctionReturn(0);
525 }
526 EXTERN_C_END
527 
528 typedef int (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
529 EXTERN_C_BEGIN
530 #undef __FUNCT__
531 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
532 int TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
533 {
534   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
535 
536   PetscFunctionBegin;
537   pseudo->dt      = dt;
538   pseudo->dtctx   = ctx;
539   PetscFunctionReturn(0);
540 }
541 EXTERN_C_END
542 
543 /* ----------------------------------------------------------------------------- */
544 
545 EXTERN_C_BEGIN
546 #undef __FUNCT__
547 #define __FUNCT__ "TSCreate_Pseudo"
548 int TSCreate_Pseudo(TS ts)
549 {
550   TS_Pseudo  *pseudo;
551   int        ierr;
552 
553   PetscFunctionBegin;
554   ts->ops->destroy         = TSDestroy_Pseudo;
555   ts->ops->view            = TSView_Pseudo;
556 
557   if (ts->problem_type == TS_LINEAR) {
558     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
559   }
560   if (!ts->A) {
561     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
562   }
563 
564   ts->ops->setup           = TSSetUp_Pseudo;
565   ts->ops->step            = TSStep_Pseudo;
566   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
567 
568   /* create the required nonlinear solver context */
569   ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr);
570 
571   ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr);
572   PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
573 
574   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
575   ts->data = (void*)pseudo;
576 
577   pseudo->dt_increment                 = 1.1;
578   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
579   pseudo->dt                           = TSPseudoDefaultTimeStep;
580 
581   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
582                     "TSPseudoSetVerifyTimeStep_Pseudo",
583                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
584   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
585                     "TSPseudoSetTimeStepIncrement_Pseudo",
586                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
587   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
588                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
589                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
590   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
591                     "TSPseudoSetTimeStep_Pseudo",
592                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 EXTERN_C_END
596 
597 #undef __FUNCT__
598 #define __FUNCT__ "TSPseudoDefaultTimeStep"
599 /*@C
600    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
601    Use with TSPseudoSetTimeStep().
602 
603    Collective on TS
604 
605    Input Parameters:
606 .  ts - the timestep context
607 .  dtctx - unused timestep context
608 
609    Output Parameter:
610 .  newdt - the timestep to use for the next step
611 
612    Level: advanced
613 
614 .keywords: timestep, pseudo, default
615 
616 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
617 @*/
618 int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
619 {
620   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
621   PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
622   int       ierr;
623 
624   PetscFunctionBegin;
625   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
626   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
627   if (pseudo->initial_fnorm == 0.0) {
628     /* first time through so compute initial function norm */
629     pseudo->initial_fnorm = pseudo->fnorm;
630     fnorm_previous        = pseudo->fnorm;
631   }
632   if (pseudo->fnorm == 0.0) {
633     *newdt = 1.e12*inc*ts->time_step;
634   } else if (pseudo->increment_dt_from_initial_dt) {
635     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
636   } else {
637     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
638   }
639   pseudo->fnorm_previous = pseudo->fnorm;
640   PetscFunctionReturn(0);
641 }
642 
643 
644 
645 
646 
647 
648 
649 
650 
651 
652 
653 
654 
655 
656 
657 
658 
659