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