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