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