xref: /petsc/src/ts/impls/pseudo/posindep.c (revision 9f95472905e43bbd86a13e9b60bb1eb5a5ab8b35)
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);CHKERRQ(ierr);
160       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
161       ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
162       ts->nonlinear_its += its; ts->linear_its += lits;
163       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
164       if (ok) break;
165       ts->ptime        -= current_time_step;
166       current_time_step = ts->time_step;
167     }
168     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
169     ts->steps++;
170   }
171   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
172   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
173   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
174 
175   *steps += ts->steps;
176   *ptime  = ts->ptime;
177   PetscFunctionReturn(0);
178 }
179 
180 /*------------------------------------------------------------*/
181 #undef __FUNCT__
182 #define __FUNCT__ "TSDestroy_Pseudo"
183 static int TSDestroy_Pseudo(TS ts)
184 {
185   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
186   int       ierr;
187 
188   PetscFunctionBegin;
189   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
190   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
191   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
192   ierr = PetscFree(pseudo);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 
197 /*------------------------------------------------------------*/
198 
199 /*
200     This defines the nonlinear equation that is to be solved with SNES
201 
202               (U^{n+1} - U^{n})/dt - F(U^{n+1})
203 */
204 #undef __FUNCT__
205 #define __FUNCT__ "TSPseudoFunction"
206 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
207 {
208   TS     ts = (TS) ctx;
209   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
210   int    ierr,i,n;
211 
212   PetscFunctionBegin;
213   /* apply user provided function */
214   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
215   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
216   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
217   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
218   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
219   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
220   for (i=0; i<n; i++) {
221     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
222   }
223   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
224   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
225   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
226 
227   PetscFunctionReturn(0);
228 }
229 
230 /*
231    This constructs the Jacobian needed for SNES
232 
233              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
234 */
235 #undef __FUNCT__
236 #define __FUNCT__ "TSPseudoJacobian"
237 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
238 {
239   TS            ts = (TS) ctx;
240   int           ierr;
241   PetscScalar   mone = -1.0,mdt = 1.0/ts->time_step;
242 
243   PetscFunctionBegin;
244   /* construct users Jacobian */
245   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
246 
247   /* shift and scale Jacobian */
248   ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
249   ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
250   if (*BB != *AA && *str != SAME_PRECONDITIONER) {
251     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
252     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
253   }
254 
255   PetscFunctionReturn(0);
256 }
257 
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "TSSetUp_Pseudo"
261 static int TSSetUp_Pseudo(TS ts)
262 {
263   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
264   int       ierr;
265 
266   PetscFunctionBegin;
267   /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */
268   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
269   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
270   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
271   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 /*------------------------------------------------------------*/
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "TSPseudoDefaultMonitor"
278 int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
279 {
280   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
281   int       ierr;
282 
283   PetscFunctionBegin;
284   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "TSSetFromOptions_Pseudo"
290 static int TSSetFromOptions_Pseudo(TS ts)
291 {
292   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
293   int        ierr;
294   PetscTruth flg;
295 
296   PetscFunctionBegin;
297 
298   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
299     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
300     if (flg) {
301       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
302     }
303     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
304     if (flg) {
305       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
306     }
307     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
308   ierr = PetscOptionsTail();CHKERRQ(ierr);
309 
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "TSView_Pseudo"
315 static int TSView_Pseudo(TS ts,PetscViewer viewer)
316 {
317   PetscFunctionBegin;
318   PetscFunctionReturn(0);
319 }
320 
321 /* ----------------------------------------------------------------------------- */
322 #undef __FUNCT__
323 #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
324 /*@C
325    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
326    last timestep.
327 
328    Collective on TS
329 
330    Input Parameters:
331 +  ts - timestep context
332 .  dt - user-defined function to verify timestep
333 -  ctx - [optional] user-defined context for private data
334          for the timestep verification routine (may be PETSC_NULL)
335 
336    Level: advanced
337 
338    Calling sequence of func:
339 .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
340 
341 .  update - latest solution vector
342 .  ctx - [optional] timestep context
343 .  newdt - the timestep to use for the next step
344 .  flag - flag indicating whether the last time step was acceptable
345 
346    Notes:
347    The routine set here will be called by TSPseudoVerifyTimeStep()
348    during the timestepping process.
349 
350 .keywords: timestep, pseudo, set, verify
351 
352 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
353 @*/
354 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
355 {
356   int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *);
357 
358   PetscFunctionBegin;
359   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
360 
361   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
362   if (f) {
363     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
364   }
365   PetscFunctionReturn(0);
366 }
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
370 /*@
371     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
372     dt when using the TSPseudoDefaultTimeStep() routine.
373 
374    Collective on TS
375 
376     Input Parameters:
377 +   ts - the timestep context
378 -   inc - the scaling factor >= 1.0
379 
380     Options Database Key:
381 $    -ts_pseudo_increment <increment>
382 
383     Level: advanced
384 
385 .keywords: timestep, pseudo, set, increment
386 
387 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
388 @*/
389 int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
390 {
391   int ierr,(*f)(TS,PetscReal);
392 
393   PetscFunctionBegin;
394   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
395 
396   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
397   if (f) {
398     ierr = (*f)(ts,inc);CHKERRQ(ierr);
399   }
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
405 /*@
406     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
407     is computed via the formula
408 $         dt = initial_dt*initial_fnorm/current_fnorm
409       rather than the default update,
410 $         dt = current_dt*previous_fnorm/current_fnorm.
411 
412    Collective on TS
413 
414     Input Parameter:
415 .   ts - the timestep context
416 
417     Options Database Key:
418 $    -ts_pseudo_increment_dt_from_initial_dt
419 
420     Level: advanced
421 
422 .keywords: timestep, pseudo, set, increment
423 
424 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
425 @*/
426 int TSPseudoIncrementDtFromInitialDt(TS ts)
427 {
428   int ierr,(*f)(TS);
429 
430   PetscFunctionBegin;
431   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
432 
433   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
434   if (f) {
435     ierr = (*f)(ts);CHKERRQ(ierr);
436   }
437   PetscFunctionReturn(0);
438 }
439 
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "TSPseudoSetTimeStep"
443 /*@C
444    TSPseudoSetTimeStep - Sets the user-defined routine to be
445    called at each pseudo-timestep to update the timestep.
446 
447    Collective on TS
448 
449    Input Parameters:
450 +  ts - timestep context
451 .  dt - function to compute timestep
452 -  ctx - [optional] user-defined context for private data
453          required by the function (may be PETSC_NULL)
454 
455    Level: intermediate
456 
457    Calling sequence of func:
458 .  func (TS ts,PetscReal *newdt,void *ctx);
459 
460 .  newdt - the newly computed timestep
461 .  ctx - [optional] timestep context
462 
463    Notes:
464    The routine set here will be called by TSPseudoComputeTimeStep()
465    during the timestepping process.
466 
467 .keywords: timestep, pseudo, set
468 
469 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
470 @*/
471 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
472 {
473   int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *);
474 
475   PetscFunctionBegin;
476   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
477 
478   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
479   if (f) {
480     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
481   }
482   PetscFunctionReturn(0);
483 }
484 
485 /* ----------------------------------------------------------------------------- */
486 
487 typedef int (*FCN1)(TS,Vec,void*,PetscReal*,int*); /* force argument to next function to not be extern C*/
488 EXTERN_C_BEGIN
489 #undef __FUNCT__
490 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
491 int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx)
492 {
493   TS_Pseudo *pseudo;
494 
495   PetscFunctionBegin;
496   pseudo              = (TS_Pseudo*)ts->data;
497   pseudo->verify      = dt;
498   pseudo->verifyctx   = ctx;
499   PetscFunctionReturn(0);
500 }
501 EXTERN_C_END
502 
503 EXTERN_C_BEGIN
504 #undef __FUNCT__
505 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo"
506 int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
507 {
508   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
509 
510   PetscFunctionBegin;
511   pseudo->dt_increment = inc;
512   PetscFunctionReturn(0);
513 }
514 EXTERN_C_END
515 
516 EXTERN_C_BEGIN
517 #undef __FUNCT__
518 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
519 int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
520 {
521   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
522 
523   PetscFunctionBegin;
524   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
525   PetscFunctionReturn(0);
526 }
527 EXTERN_C_END
528 
529 typedef int (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
530 EXTERN_C_BEGIN
531 #undef __FUNCT__
532 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
533 int TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
534 {
535   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
536 
537   PetscFunctionBegin;
538   pseudo->dt      = dt;
539   pseudo->dtctx   = ctx;
540   PetscFunctionReturn(0);
541 }
542 EXTERN_C_END
543 
544 /* ----------------------------------------------------------------------------- */
545 
546 EXTERN_C_BEGIN
547 #undef __FUNCT__
548 #define __FUNCT__ "TSCreate_Pseudo"
549 int TSCreate_Pseudo(TS ts)
550 {
551   TS_Pseudo  *pseudo;
552   int        ierr;
553 
554   PetscFunctionBegin;
555   ts->ops->destroy         = TSDestroy_Pseudo;
556   ts->ops->view            = TSView_Pseudo;
557 
558   if (ts->problem_type == TS_LINEAR) {
559     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
560   }
561   if (!ts->A) {
562     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian");
563   }
564 
565   ts->ops->setup           = TSSetUp_Pseudo;
566   ts->ops->step            = TSStep_Pseudo;
567   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
568 
569   /* create the required nonlinear solver context */
570   ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr);
571 
572   ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr);
573   PetscLogObjectMemory(ts,sizeof(TS_Pseudo));
574 
575   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
576   ts->data = (void*)pseudo;
577 
578   pseudo->dt_increment                 = 1.1;
579   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
580   pseudo->dt                           = TSPseudoDefaultTimeStep;
581 
582   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
583                     "TSPseudoSetVerifyTimeStep_Pseudo",
584                      TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
585   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
586                     "TSPseudoSetTimeStepIncrement_Pseudo",
587                      TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
588   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
589                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
590                      TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
591   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
592                     "TSPseudoSetTimeStep_Pseudo",
593                      TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
594   PetscFunctionReturn(0);
595 }
596 EXTERN_C_END
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "TSPseudoDefaultTimeStep"
600 /*@C
601    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
602    Use with TSPseudoSetTimeStep().
603 
604    Collective on TS
605 
606    Input Parameters:
607 .  ts - the timestep context
608 .  dtctx - unused timestep context
609 
610    Output Parameter:
611 .  newdt - the timestep to use for the next step
612 
613    Level: advanced
614 
615 .keywords: timestep, pseudo, default
616 
617 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
618 @*/
619 int TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
620 {
621   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
622   PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
623   int       ierr;
624 
625   PetscFunctionBegin;
626   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
627   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
628   if (pseudo->initial_fnorm == 0.0) {
629     /* first time through so compute initial function norm */
630     pseudo->initial_fnorm = pseudo->fnorm;
631     fnorm_previous        = pseudo->fnorm;
632   }
633   if (pseudo->fnorm == 0.0) {
634     *newdt = 1.e12*inc*ts->time_step;
635   } else if (pseudo->increment_dt_from_initial_dt) {
636     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
637   } else {
638     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
639   }
640   pseudo->fnorm_previous = pseudo->fnorm;
641   PetscFunctionReturn(0);
642 }
643 
644 
645 
646 
647 
648 
649 
650 
651 
652 
653 
654 
655 
656 
657 
658 
659 
660