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