xref: /petsc/src/ts/impls/pseudo/posindep.c (revision b1d4fb267e72f6b087ed013bbfbbee1418c3f503)
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 EXTERN_C_BEGIN
8 typedef struct {
9   Vec  update;      /* work vector where new solution is formed */
10   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
11   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
12 
13   /* information used for Pseudo-timestepping */
14 
15   int    (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
16   void   *dtctx;
17   int    (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */
18   void   *verifyctx;
19 
20   PetscReal  initial_fnorm,fnorm;                  /* original and current norm of F(u) */
21   PetscReal  fnorm_previous;
22 
23   PetscReal  dt_increment;                  /* scaling that dt is incremented each time-step */
24   PetscTruth increment_dt_from_initial_dt;
25 } TS_Pseudo;
26 EXTERN_C_END
27 
28 /* ------------------------------------------------------------------------------*/
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "TSPseudoComputeTimeStep"
32 /*@
33     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
34     pseudo-timestepping process.
35 
36     Collective on TS
37 
38     Input Parameter:
39 .   ts - timestep context
40 
41     Output Parameter:
42 .   dt - newly computed timestep
43 
44     Level: advanced
45 
46     Notes:
47     The routine to be called here to compute the timestep should be
48     set by calling TSPseudoSetTimeStep().
49 
50 .keywords: timestep, pseudo, compute
51 
52 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep()
53 @*/
54 int TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
55 {
56   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
57   int       ierr;
58 
59   PetscFunctionBegin;
60   ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
61   ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr);
62   ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 
67 /* ------------------------------------------------------------------------------*/
68 #undef __FUNCT__
69 #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep"
70 /*@C
71    TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep.
72 
73    Collective on TS
74 
75    Input Parameters:
76 +  ts - the timestep context
77 .  dtctx - unused timestep context
78 -  update - latest solution vector
79 
80    Output Parameters:
81 +  newdt - the timestep to use for the next step
82 -  flag - flag indicating whether the last time step was acceptable
83 
84    Level: advanced
85 
86    Note:
87    This routine always returns a flag of 1, indicating an acceptable
88    timestep.
89 
90 .keywords: timestep, pseudo, default, verify
91 
92 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
93 @*/
94 int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag)
95 {
96   PetscFunctionBegin;
97   *flag = 1;
98   PetscFunctionReturn(0);
99 }
100 
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "TSPseudoVerifyTimeStep"
104 /*@
105     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
106 
107     Collective on TS
108 
109     Input Parameters:
110 +   ts - timestep context
111 -   update - latest solution vector
112 
113     Output Parameters:
114 +   dt - newly computed timestep (if it had to shrink)
115 -   flag - indicates if current timestep was ok
116 
117     Level: advanced
118 
119     Notes:
120     The routine to be called here to compute the timestep should be
121     set by calling TSPseudoSetVerifyTimeStep().
122 
123 .keywords: timestep, pseudo, verify
124 
125 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep()
126 @*/
127 int TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag)
128 {
129   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
130   int       ierr;
131 
132   PetscFunctionBegin;
133   if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);}
134 
135   ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr);
136 
137   PetscFunctionReturn(0);
138 }
139 
140 /* --------------------------------------------------------------------------------*/
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "TSStep_Pseudo"
144 static int TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime)
145 {
146   Vec       sol = ts->vec_sol;
147   int       ierr,i,max_steps = ts->max_steps,its,ok,lits;
148   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
149   PetscReal current_time_step;
150 
151   PetscFunctionBegin;
152   *steps = -ts->steps;
153 
154   ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr);
155   for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) {
156     ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr);
157     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
158     current_time_step = ts->time_step;
159     while (PETSC_TRUE) {
160       ts->ptime  += current_time_step;
161       ierr = SNESSolve(ts->snes,pseudo->update,&its);CHKERRQ(ierr);
162       ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr);
163       ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits;
164       ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr);
165       if (ok) break;
166       ts->ptime        -= current_time_step;
167       current_time_step = ts->time_step;
168     }
169     ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr);
170     ts->steps++;
171   }
172   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
173   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
174   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
175 
176   *steps += ts->steps;
177   *ptime  = ts->ptime;
178   PetscFunctionReturn(0);
179 }
180 
181 /*------------------------------------------------------------*/
182 #undef __FUNCT__
183 #define __FUNCT__ "TSDestroy_Pseudo"
184 static int TSDestroy_Pseudo(TS ts)
185 {
186   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
187   int       ierr;
188 
189   PetscFunctionBegin;
190   if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);}
191   if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);}
192   if (pseudo->rhs)  {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);}
193   ierr = PetscFree(pseudo);CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 
198 /*------------------------------------------------------------*/
199 
200 /*
201     This defines the nonlinear equation that is to be solved with SNES
202 
203               (U^{n+1} - U^{n})/dt - F(U^{n+1})
204 */
205 #undef __FUNCT__
206 #define __FUNCT__ "TSPseudoFunction"
207 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx)
208 {
209   TS     ts = (TS) ctx;
210   PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
211   int    ierr,i,n;
212 
213   PetscFunctionBegin;
214   /* apply user provided function */
215   ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr);
216   /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */
217   ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr);
218   ierr = VecGetArray(x,&unp1);CHKERRQ(ierr);
219   ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr);
220   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
221   for (i=0; i<n; i++) {
222     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
223   }
224   ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr);
225   ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr);
226   ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr);
227 
228   PetscFunctionReturn(0);
229 }
230 
231 /*
232    This constructs the Jacobian needed for SNES
233 
234              J = I/dt - J_{F}   where J_{F} is the given Jacobian of F.
235 */
236 #undef __FUNCT__
237 #define __FUNCT__ "TSPseudoJacobian"
238 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
239 {
240   TS            ts = (TS) ctx;
241   int           ierr;
242   PetscScalar   mone = -1.0,mdt = 1.0/ts->time_step;
243 
244   PetscFunctionBegin;
245   /* construct users Jacobian */
246   ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr);
247 
248   /* shift and scale Jacobian */
249   ierr = MatScale(&mone,*AA);CHKERRQ(ierr);
250   ierr = MatShift(&mdt,*AA);CHKERRQ(ierr);
251   if (*BB != *AA && *str != SAME_PRECONDITIONER) {
252     ierr = MatScale(&mone,*BB);CHKERRQ(ierr);
253     ierr = MatShift(&mdt,*BB);CHKERRQ(ierr);
254   }
255 
256   PetscFunctionReturn(0);
257 }
258 
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "TSSetUp_Pseudo"
262 static int TSSetUp_Pseudo(TS ts)
263 {
264   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
265   int       ierr;
266 
267   PetscFunctionBegin;
268   /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */
269   ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr);
270   ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr);
271   ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr);
272   ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr);
273   PetscFunctionReturn(0);
274 }
275 /*------------------------------------------------------------*/
276 
277 #undef __FUNCT__
278 #define __FUNCT__ "TSPseudoDefaultMonitor"
279 int TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
280 {
281   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
282   int       ierr;
283 
284   PetscFunctionBegin;
285   ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "TSSetFromOptions_Pseudo"
291 static int TSSetFromOptions_Pseudo(TS ts)
292 {
293   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
294   int        ierr;
295   PetscTruth flg;
296 
297   PetscFunctionBegin;
298 
299   ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr);
300     ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr);
301     if (flg) {
302       ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
303     }
304     ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr);
305     if (flg) {
306       ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
307     }
308     ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr);
309   ierr = PetscOptionsTail();CHKERRQ(ierr);
310 
311   PetscFunctionReturn(0);
312 }
313 
314 #undef __FUNCT__
315 #define __FUNCT__ "TSView_Pseudo"
316 static int TSView_Pseudo(TS ts,PetscViewer viewer)
317 {
318   PetscFunctionBegin;
319   PetscFunctionReturn(0);
320 }
321 
322 /* ----------------------------------------------------------------------------- */
323 #undef __FUNCT__
324 #define __FUNCT__ "TSPseudoSetVerifyTimeStep"
325 /*@C
326    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
327    last timestep.
328 
329    Collective on TS
330 
331    Input Parameters:
332 +  ts - timestep context
333 .  dt - user-defined function to verify timestep
334 -  ctx - [optional] user-defined context for private data
335          for the timestep verification routine (may be PETSC_NULL)
336 
337    Level: advanced
338 
339    Calling sequence of func:
340 .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag);
341 
342 .  update - latest solution vector
343 .  ctx - [optional] timestep context
344 .  newdt - the timestep to use for the next step
345 .  flag - flag indicating whether the last time step was acceptable
346 
347    Notes:
348    The routine set here will be called by TSPseudoVerifyTimeStep()
349    during the timestepping process.
350 
351 .keywords: timestep, pseudo, set, verify
352 
353 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
354 @*/
355 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx)
356 {
357   int ierr,(*f)(TS,int (*)(TS,Vec,void*,PetscReal *,int *),void *);
358 
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(ts,TS_COOKIE);
361 
362   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
363   if (f) {
364     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
365   }
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "TSPseudoSetTimeStepIncrement"
371 /*@
372     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
373     dt when using the TSPseudoDefaultTimeStep() routine.
374 
375    Collective on TS
376 
377     Input Parameters:
378 +   ts - the timestep context
379 -   inc - the scaling factor >= 1.0
380 
381     Options Database Key:
382 $    -ts_pseudo_increment <increment>
383 
384     Level: advanced
385 
386 .keywords: timestep, pseudo, set, increment
387 
388 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
389 @*/
390 int TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
391 {
392   int ierr,(*f)(TS,PetscReal);
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(ts,TS_COOKIE);
396 
397   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr);
398   if (f) {
399     ierr = (*f)(ts,inc);CHKERRQ(ierr);
400   }
401   PetscFunctionReturn(0);
402 }
403 
404 #undef __FUNCT__
405 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt"
406 /*@
407     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
408     is computed via the formula
409 $         dt = initial_dt*initial_fnorm/current_fnorm
410       rather than the default update,
411 $         dt = current_dt*previous_fnorm/current_fnorm.
412 
413    Collective on TS
414 
415     Input Parameter:
416 .   ts - the timestep context
417 
418     Options Database Key:
419 $    -ts_pseudo_increment_dt_from_initial_dt
420 
421     Level: advanced
422 
423 .keywords: timestep, pseudo, set, increment
424 
425 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
426 @*/
427 int TSPseudoIncrementDtFromInitialDt(TS ts)
428 {
429   int ierr,(*f)(TS);
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(ts,TS_COOKIE);
433 
434   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr);
435   if (f) {
436     ierr = (*f)(ts);CHKERRQ(ierr);
437   }
438   PetscFunctionReturn(0);
439 }
440 
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "TSPseudoSetTimeStep"
444 /*@C
445    TSPseudoSetTimeStep - Sets the user-defined routine to be
446    called at each pseudo-timestep to update the timestep.
447 
448    Collective on TS
449 
450    Input Parameters:
451 +  ts - timestep context
452 .  dt - function to compute timestep
453 -  ctx - [optional] user-defined context for private data
454          required by the function (may be PETSC_NULL)
455 
456    Level: intermediate
457 
458    Calling sequence of func:
459 .  func (TS ts,PetscReal *newdt,void *ctx);
460 
461 .  newdt - the newly computed timestep
462 .  ctx - [optional] timestep context
463 
464    Notes:
465    The routine set here will be called by TSPseudoComputeTimeStep()
466    during the timestepping process.
467 
468 .keywords: timestep, pseudo, set
469 
470 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
471 @*/
472 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,PetscReal*,void*),void* ctx)
473 {
474   int ierr,(*f)(TS,int (*)(TS,PetscReal *,void *),void *);
475 
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(ts,TS_COOKIE);
478 
479   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr);
480   if (f) {
481     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
482   }
483   PetscFunctionReturn(0);
484 }
485 
486 /* ----------------------------------------------------------------------------- */
487 
488 EXTERN_C_BEGIN
489 #undef __FUNCT__
490 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo"
491 int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,PetscReal*,int*),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 EXTERN_C_BEGIN
530 #undef __FUNCT__
531 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo"
532 int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,PetscReal*,void*),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