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