xref: /petsc/src/ts/impls/pseudo/posindep.c (revision f1af5d2ffeae1f5fc391a89939f4818e47770ae3)
1 /*$Id: posindep.c,v 1.37 1999/10/24 14:03:53 bsmith 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;
324   PetscTruth flg;
325   double     inc;
326 
327   PetscFunctionBegin;
328   ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
329 
330   ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg);CHKERRQ(ierr);
331   if (flg) {
332     ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0);CHKERRQ(ierr);
333   }
334   ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg);CHKERRQ(ierr);
335   if (flg) {
336     ierr = TSPseudoSetTimeStepIncrement(ts,inc);CHKERRQ(ierr);
337   }
338   ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr);
339   if (flg) {
340     ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr);
341   }
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNC__
346 #define __FUNC__ "TSPrintHelp_Pseudo"
347 static int TSPrintHelp_Pseudo(TS ts,char *p)
348 {
349   int ierr;
350 
351   PetscFunctionBegin;
352   ierr = (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n");CHKERRQ(ierr);
353   ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p);CHKERRQ(ierr);
354   ierr = (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p);CHKERRQ(ierr);
355   ierr = (*PetscHelpPrintf)(ts->comm,"     initial fnorm/current fnorm to determine new timestep\n");CHKERRQ(ierr);
356   ierr = (*PetscHelpPrintf)(ts->comm,"     default is current_dt * previous fnorm/current fnorm\n");CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNC__
361 #define __FUNC__ "TSView_Pseudo"
362 static int TSView_Pseudo(TS ts,Viewer viewer)
363 {
364   PetscFunctionBegin;
365   PetscFunctionReturn(0);
366 }
367 
368 /* ----------------------------------------------------------------------------- */
369 #undef __FUNC__
370 #define __FUNC__ "TSPseudoSetVerifyTimeStep"
371 /*@
372    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
373    last timestep.
374 
375    Collective on TS
376 
377    Input Parameters:
378 +  ts - timestep context
379 .  dt - user-defined function to verify timestep
380 -  ctx - [optional] user-defined context for private data
381          for the timestep verification routine (may be PETSC_NULL)
382 
383    Level: advanced
384 
385    Calling sequence of func:
386 .  func (TS ts,Vec update,void *ctx,double *newdt,int *flag);
387 
388 .  update - latest solution vector
389 .  ctx - [optional] timestep context
390 .  newdt - the timestep to use for the next step
391 .  flag - flag indicating whether the last time step was acceptable
392 
393    Notes:
394    The routine set here will be called by TSPseudoVerifyTimeStep()
395    during the timestepping process.
396 
397 .keywords: timestep, pseudo, set, verify
398 
399 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep()
400 @*/
401 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
402 {
403   int ierr, (*f)(TS,int (*)(TS,Vec,void*,double *,int *),void *);
404 
405   PetscFunctionBegin;
406   PetscValidHeaderSpecific(ts,TS_COOKIE);
407 
408   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void **)&f);CHKERRQ(ierr);
409   if (f) {
410     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
411   }
412   PetscFunctionReturn(0);
413 }
414 
415 #undef __FUNC__
416 #define __FUNC__ "TSPseudoSetTimeStepIncrement"
417 /*@
418     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
419     dt when using the TSPseudoDefaultTimeStep() routine.
420 
421    Collective on TS
422 
423     Input Parameters:
424 +   ts - the timestep context
425 -   inc - the scaling factor >= 1.0
426 
427     Options Database Key:
428 $    -ts_pseudo_increment <increment>
429 
430     Level: advanced
431 
432 .keywords: timestep, pseudo, set, increment
433 
434 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
435 @*/
436 int TSPseudoSetTimeStepIncrement(TS ts,double inc)
437 {
438   int ierr, (*f)(TS,double);
439 
440   PetscFunctionBegin;
441   PetscValidHeaderSpecific(ts,TS_COOKIE);
442 
443   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void **)&f);CHKERRQ(ierr);
444   if (f) {
445     ierr = (*f)(ts,inc);CHKERRQ(ierr);
446   }
447   PetscFunctionReturn(0);
448 }
449 
450 #undef __FUNC__
451 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt"
452 /*@
453     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
454     is computed via the formula
455 $         dt = initial_dt*initial_fnorm/current_fnorm
456       rather than the default update,
457 $         dt = current_dt*previous_fnorm/current_fnorm.
458 
459    Collective on TS
460 
461     Input Parameter:
462 .   ts - the timestep context
463 
464     Options Database Key:
465 $    -ts_pseudo_increment_dt_from_initial_dt
466 
467     Level: advanced
468 
469 .keywords: timestep, pseudo, set, increment
470 
471 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep()
472 @*/
473 int TSPseudoIncrementDtFromInitialDt(TS ts)
474 {
475   int ierr, (*f)(TS);
476 
477   PetscFunctionBegin;
478   PetscValidHeaderSpecific(ts,TS_COOKIE);
479 
480   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void **)&f);CHKERRQ(ierr);
481   if (f) {
482     ierr = (*f)(ts);CHKERRQ(ierr);
483   }
484   PetscFunctionReturn(0);
485 }
486 
487 
488 #undef __FUNC__
489 #define __FUNC__ "TSPseudoSetTimeStep"
490 /*@
491    TSPseudoSetTimeStep - Sets the user-defined routine to be
492    called at each pseudo-timestep to update the timestep.
493 
494    Collective on TS
495 
496    Input Parameters:
497 +  ts - timestep context
498 .  dt - function to compute timestep
499 -  ctx - [optional] user-defined context for private data
500          required by the function (may be PETSC_NULL)
501 
502    Level: intermediate
503 
504    Calling sequence of func:
505 .  func (TS ts,double *newdt,void *ctx);
506 
507 .  newdt - the newly computed timestep
508 .  ctx - [optional] timestep context
509 
510    Notes:
511    The routine set here will be called by TSPseudoComputeTimeStep()
512    during the timestepping process.
513 
514 .keywords: timestep, pseudo, set
515 
516 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
517 @*/
518 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx)
519 {
520   int ierr, (*f)(TS,int (*)(TS,double *,void *),void *);
521 
522   PetscFunctionBegin;
523   PetscValidHeaderSpecific(ts,TS_COOKIE);
524 
525   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void **)&f);CHKERRQ(ierr);
526   if (f) {
527     ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr);
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 /* ----------------------------------------------------------------------------- */
533 
534 EXTERN_C_BEGIN
535 #undef __FUNC__
536 #define __FUNC__ "TSPseudoSetVerifyTimeStep_Pseudo"
537 int TSPseudoSetVerifyTimeStep_Pseudo(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx)
538 {
539   TS_Pseudo *pseudo;
540 
541   PetscFunctionBegin;
542   pseudo              = (TS_Pseudo*) ts->data;
543   pseudo->verify      = dt;
544   pseudo->verifyctx   = ctx;
545   PetscFunctionReturn(0);
546 }
547 EXTERN_C_END
548 
549 EXTERN_C_BEGIN
550 #undef __FUNC__
551 #define __FUNC__ "TSPseudoSetTimeStepIncrement_Pseudo"
552 int TSPseudoSetTimeStepIncrement_Pseudo(TS ts,double inc)
553 {
554   TS_Pseudo *pseudo;
555 
556   PetscFunctionBegin;
557   pseudo               = (TS_Pseudo*) ts->data;
558   pseudo->dt_increment = inc;
559   PetscFunctionReturn(0);
560 }
561 EXTERN_C_END
562 
563 EXTERN_C_BEGIN
564 #undef __FUNC__
565 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt_Pseudo"
566 int TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
567 {
568   TS_Pseudo *pseudo;
569 
570   PetscFunctionBegin;
571   pseudo                               = (TS_Pseudo*) ts->data;
572   pseudo->increment_dt_from_initial_dt = 1;
573   PetscFunctionReturn(0);
574 }
575 EXTERN_C_END
576 
577 EXTERN_C_BEGIN
578 #undef __FUNC__
579 #define __FUNC__ "TSPseudoSetTimeStep_Pseudo"
580 int TSPseudoSetTimeStep_Pseudo(TS ts,int (*dt)(TS,double*,void*),void* ctx)
581 {
582   TS_Pseudo *pseudo;
583 
584   PetscFunctionBegin;
585   pseudo          = (TS_Pseudo*) ts->data;
586   pseudo->dt      = dt;
587   pseudo->dtctx   = ctx;
588   PetscFunctionReturn(0);
589 }
590 EXTERN_C_END
591 
592 /* ----------------------------------------------------------------------------- */
593 
594 EXTERN_C_BEGIN
595 #undef __FUNC__
596 #define __FUNC__ "TSCreate_Pseudo"
597 int TSCreate_Pseudo(TS ts )
598 {
599   TS_Pseudo *pseudo;
600   int       ierr;
601   MatType   mtype;
602 
603   PetscFunctionBegin;
604   ts->destroy         = TSDestroy_Pseudo;
605   ts->printhelp       = TSPrintHelp_Pseudo;
606   ts->view            = TSView_Pseudo;
607 
608   if (ts->problem_type == TS_LINEAR) {
609     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems");
610   }
611   if (!ts->A) {
612     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian");
613   }
614   ierr = MatGetType(ts->A,&mtype,PETSC_NULL);CHKERRQ(ierr);
615   if (mtype == MATSHELL) {
616     ts->Ashell = ts->A;
617   }
618   ts->setup           = TSSetUp_Pseudo;
619   ts->step            = TSStep_Pseudo;
620   ts->setfromoptions  = TSSetFromOptions_Pseudo;
621 
622   /* create the required nonlinear solver context */
623   ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr);
624 
625   pseudo   = PetscNew(TS_Pseudo);CHKPTRQ(pseudo);
626   PLogObjectMemory(ts,sizeof(TS_Pseudo));
627 
628   ierr     = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr);
629   ts->data = (void *) pseudo;
630 
631   pseudo->dt_increment                 = 1.1;
632   pseudo->increment_dt_from_initial_dt = 0;
633   pseudo->dt                           = TSPseudoDefaultTimeStep;
634 
635   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
636                     "TSPseudoSetVerifyTimeStep_Pseudo",
637                     (void*)TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr);
638   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
639                     "TSPseudoSetTimeStepIncrement_Pseudo",
640                     (void*)TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr);
641   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
642                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
643                     (void*)TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr);
644   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","TSPseudoSetTimeStep_Pseudo",
645                     (void*)TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr);
646   PetscFunctionReturn(0);
647 }
648 EXTERN_C_END
649 
650 #undef __FUNC__
651 #define __FUNC__ "TSPseudoDefaultTimeStep"
652 /*@C
653    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
654    Use with TSPseudoSetTimeStep().
655 
656    Collective on TS
657 
658    Input Parameters:
659 .  ts - the timestep context
660 .  dtctx - unused timestep context
661 
662    Output Parameter:
663 .  newdt - the timestep to use for the next step
664 
665    Level: advanced
666 
667 .keywords: timestep, pseudo, default
668 
669 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
670 @*/
671 int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx)
672 {
673   TS_Pseudo *pseudo = (TS_Pseudo*) ts->data;
674   double    inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
675   int       ierr;
676 
677   PetscFunctionBegin;
678   ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr);
679   ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr);
680   if (pseudo->initial_fnorm == 0.0) {
681     /* first time through so compute initial function norm */
682     pseudo->initial_fnorm = pseudo->fnorm;
683     fnorm_previous        = pseudo->fnorm;
684   }
685   if (pseudo->fnorm == 0.0) {
686     *newdt = 1.e12*inc*ts->time_step;
687   }
688   else if (pseudo->increment_dt_from_initial_dt) {
689     *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm;
690   } else {
691     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
692   }
693   pseudo->fnorm_previous = pseudo->fnorm;
694   PetscFunctionReturn(0);
695 }
696 
697