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