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