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