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