Actual source code: posindep.c
petsc-3.6.4 2016-04-12
1: /*
2: Code for Timestepping with implicit backwards Euler.
3: */
4: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
6: typedef struct {
7: Vec update; /* work vector where new solution is formed */
8: Vec func; /* work vector where F(t[i],u[i]) is stored */
9: Vec xdot; /* work vector for time derivative of state */
11: /* information used for Pseudo-timestepping */
13: PetscErrorCode (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */
14: void *dtctx;
15: PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscBool*); /* verify previous timestep and related context */
16: void *verifyctx;
18: PetscReal fnorm_initial,fnorm; /* original and current norm of F(u) */
19: PetscReal fnorm_previous;
21: PetscReal dt_initial; /* initial time-step */
22: PetscReal dt_increment; /* scaling that dt is incremented each time-step */
23: PetscReal dt_max; /* maximum time step */
24: PetscBool increment_dt_from_initial_dt;
25: } TS_Pseudo;
27: /* ------------------------------------------------------------------------------*/
31: /*@C
32: TSPseudoComputeTimeStep - Computes the next timestep for a currently running
33: pseudo-timestepping process.
35: Collective on TS
37: Input Parameter:
38: . ts - timestep context
40: Output Parameter:
41: . dt - newly computed timestep
43: Level: developer
45: Notes:
46: The routine to be called here to compute the timestep should be
47: set by calling TSPseudoSetTimeStep().
49: .keywords: timestep, pseudo, compute
51: .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep()
52: @*/
53: PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
54: {
55: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
59: PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
60: (*pseudo->dt)(ts,dt,pseudo->dtctx);
61: PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
62: return(0);
63: }
66: /* ------------------------------------------------------------------------------*/
69: /*@C
70: TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.
72: Collective on TS
74: Input Parameters:
75: + ts - the timestep context
76: . dtctx - unused timestep context
77: - update - latest solution vector
79: Output Parameters:
80: + newdt - the timestep to use for the next step
81: - flag - flag indicating whether the last time step was acceptable
83: Level: advanced
85: Note:
86: This routine always returns a flag of 1, indicating an acceptable
87: timestep.
89: .keywords: timestep, pseudo, default, verify
91: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
92: @*/
93: PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag)
94: {
96: *flag = PETSC_TRUE;
97: return(0);
98: }
103: /*@
104: TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.
106: Collective on TS
108: Input Parameters:
109: + ts - timestep context
110: - update - latest solution vector
112: Output Parameters:
113: + dt - newly computed timestep (if it had to shrink)
114: - flag - indicates if current timestep was ok
116: Level: advanced
118: Notes:
119: The routine to be called here to compute the timestep should be
120: set by calling TSPseudoSetVerifyTimeStep().
122: .keywords: timestep, pseudo, verify
124: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
125: @*/
126: PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag)
127: {
128: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
132: if (!pseudo->verify) {*flag = PETSC_TRUE; return(0);}
134: (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);
135: return(0);
136: }
138: /* --------------------------------------------------------------------------------*/
142: static PetscErrorCode TSStep_Pseudo(TS ts)
143: {
144: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
145: PetscInt its,lits,reject;
146: PetscBool stepok;
147: PetscReal next_time_step;
148: SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
149: PetscErrorCode ierr;
152: if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
153: VecCopy(ts->vec_sol,pseudo->update);
154: next_time_step = ts->time_step;
155: TSPseudoComputeTimeStep(ts,&next_time_step);
156: for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
157: ts->time_step = next_time_step;
158: TSPreStep(ts);
159: TSPreStage(ts,ts->ptime+ts->time_step);
160: SNESSolve(ts->snes,NULL,pseudo->update);
161: SNESGetConvergedReason(ts->snes,&snesreason);
162: SNESGetLinearSolveIterations(ts->snes,&lits);
163: SNESGetIterationNumber(ts->snes,&its);
164: TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));
165: ts->snes_its += its; ts->ksp_its += lits;
166: PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);
167: pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */
168: TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);
169: if (stepok) break;
170: }
171: if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
172: ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
173: PetscInfo2(ts,"step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
174: return(0);
175: }
176: if (reject >= ts->max_reject) {
177: ts->reason = TS_DIVERGED_STEP_REJECTED;
178: PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
179: return(0);
180: }
181: VecCopy(pseudo->update,ts->vec_sol);
182: ts->ptime += ts->time_step;
183: ts->time_step = next_time_step;
184: ts->steps++;
185: return(0);
186: }
188: /*------------------------------------------------------------*/
191: static PetscErrorCode TSReset_Pseudo(TS ts)
192: {
193: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
197: VecDestroy(&pseudo->update);
198: VecDestroy(&pseudo->func);
199: VecDestroy(&pseudo->xdot);
200: return(0);
201: }
205: static PetscErrorCode TSDestroy_Pseudo(TS ts)
206: {
210: TSReset_Pseudo(ts);
211: PetscFree(ts->data);
212: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);
213: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);
214: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);
215: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);
216: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);
217: return(0);
218: }
220: /*------------------------------------------------------------*/
224: /*
225: Compute Xdot = (X^{n+1}-X^n)/dt) = 0
226: */
227: static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
228: {
229: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
230: const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn;
231: PetscScalar *xdot;
232: PetscErrorCode ierr;
233: PetscInt i,n;
236: VecGetArrayRead(ts->vec_sol,&xn);
237: VecGetArrayRead(X,&xnp1);
238: VecGetArray(pseudo->xdot,&xdot);
239: VecGetLocalSize(X,&n);
240: for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
241: VecRestoreArrayRead(ts->vec_sol,&xn);
242: VecRestoreArrayRead(X,&xnp1);
243: VecRestoreArray(pseudo->xdot,&xdot);
244: *Xdot = pseudo->xdot;
245: return(0);
246: }
250: /*
251: The transient residual is
253: F(U^{n+1},(U^{n+1}-U^n)/dt) = 0
255: or for ODE,
257: (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0
259: This is the function that must be evaluated for transient simulation and for
260: finite difference Jacobians. On the first Newton step, this algorithm uses
261: a guess of U^{n+1} = U^n in which case the transient term vanishes and the
262: residual is actually the steady state residual. Pseudotransient
263: continuation as described in the literature is a linearly implicit
264: algorithm, it only takes this one Newton step with the steady state
265: residual, and then advances to the next time step.
266: */
267: static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
268: {
269: Vec Xdot;
273: TSPseudoGetXdot(ts,X,&Xdot);
274: TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);
275: return(0);
276: }
280: /*
281: This constructs the Jacobian needed for SNES. For DAE, this is
283: dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot
285: and for ODE:
287: J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs.
288: */
289: static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
290: {
291: Vec Xdot;
295: TSPseudoGetXdot(ts,X,&Xdot);
296: TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);
297: return(0);
298: }
303: static PetscErrorCode TSSetUp_Pseudo(TS ts)
304: {
305: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
309: VecDuplicate(ts->vec_sol,&pseudo->update);
310: VecDuplicate(ts->vec_sol,&pseudo->func);
311: VecDuplicate(ts->vec_sol,&pseudo->xdot);
312: return(0);
313: }
314: /*------------------------------------------------------------*/
318: PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
319: {
320: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
322: PetscViewer viewer = (PetscViewer) dummy;
325: if (!viewer) {
326: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
327: }
328: if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */
329: VecZeroEntries(pseudo->xdot);
330: TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);
331: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
332: }
333: PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
334: PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);
335: PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
336: return(0);
337: }
341: static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptions *PetscOptionsObject,TS ts)
342: {
343: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
345: PetscBool flg = PETSC_FALSE;
346: PetscViewer viewer;
349: PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");
350: PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);
351: if (flg) {
352: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);
353: TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
354: }
355: flg = PETSC_FALSE;
356: PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);
357: if (flg) {
358: TSPseudoIncrementDtFromInitialDt(ts);
359: }
360: PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);
361: PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);
362: PetscOptionsTail();
363: return(0);
364: }
368: static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
369: {
373: SNESView(ts->snes,viewer);
374: return(0);
375: }
377: /* ----------------------------------------------------------------------------- */
380: /*@C
381: TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
382: last timestep.
384: Logically Collective on TS
386: Input Parameters:
387: + ts - timestep context
388: . dt - user-defined function to verify timestep
389: - ctx - [optional] user-defined context for private data
390: for the timestep verification routine (may be NULL)
392: Level: advanced
394: Calling sequence of func:
395: . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag);
397: . update - latest solution vector
398: . ctx - [optional] timestep context
399: . newdt - the timestep to use for the next step
400: . flag - flag indicating whether the last time step was acceptable
402: Notes:
403: The routine set here will be called by TSPseudoVerifyTimeStep()
404: during the timestepping process.
406: .keywords: timestep, pseudo, set, verify
408: .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
409: @*/
410: PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
411: {
416: PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));
417: return(0);
418: }
422: /*@
423: TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
424: dt when using the TSPseudoTimeStepDefault() routine.
426: Logically Collective on TS
428: Input Parameters:
429: + ts - the timestep context
430: - inc - the scaling factor >= 1.0
432: Options Database Key:
433: . -ts_pseudo_increment <increment>
435: Level: advanced
437: .keywords: timestep, pseudo, set, increment
439: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
440: @*/
441: PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
442: {
448: PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));
449: return(0);
450: }
454: /*@
455: TSPseudoSetMaxTimeStep - Sets the maximum time step
456: when using the TSPseudoTimeStepDefault() routine.
458: Logically Collective on TS
460: Input Parameters:
461: + ts - the timestep context
462: - maxdt - the maximum time step, use a non-positive value to deactivate
464: Options Database Key:
465: . -ts_pseudo_max_dt <increment>
467: Level: advanced
469: .keywords: timestep, pseudo, set
471: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
472: @*/
473: PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
474: {
480: PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
481: return(0);
482: }
486: /*@
487: TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
488: is computed via the formula
489: $ dt = initial_dt*initial_fnorm/current_fnorm
490: rather than the default update,
491: $ dt = current_dt*previous_fnorm/current_fnorm.
493: Logically Collective on TS
495: Input Parameter:
496: . ts - the timestep context
498: Options Database Key:
499: . -ts_pseudo_increment_dt_from_initial_dt
501: Level: advanced
503: .keywords: timestep, pseudo, set, increment
505: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
506: @*/
507: PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts)
508: {
513: PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));
514: return(0);
515: }
520: /*@C
521: TSPseudoSetTimeStep - Sets the user-defined routine to be
522: called at each pseudo-timestep to update the timestep.
524: Logically Collective on TS
526: Input Parameters:
527: + ts - timestep context
528: . dt - function to compute timestep
529: - ctx - [optional] user-defined context for private data
530: required by the function (may be NULL)
532: Level: intermediate
534: Calling sequence of func:
535: . func (TS ts,PetscReal *newdt,void *ctx);
537: . newdt - the newly computed timestep
538: . ctx - [optional] timestep context
540: Notes:
541: The routine set here will be called by TSPseudoComputeTimeStep()
542: during the timestepping process.
543: If not set then TSPseudoTimeStepDefault() is automatically used
545: .keywords: timestep, pseudo, set
547: .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
548: @*/
549: PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
550: {
555: PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));
556: return(0);
557: }
559: /* ----------------------------------------------------------------------------- */
561: typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/
564: PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
565: {
566: TS_Pseudo *pseudo;
569: pseudo = (TS_Pseudo*)ts->data;
570: pseudo->verify = dt;
571: pseudo->verifyctx = ctx;
572: return(0);
573: }
577: PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
578: {
579: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
582: pseudo->dt_increment = inc;
583: return(0);
584: }
588: PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
589: {
590: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
593: pseudo->dt_max = maxdt;
594: return(0);
595: }
599: PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
600: {
601: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
604: pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
605: return(0);
606: }
608: typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
611: PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
612: {
613: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
616: pseudo->dt = dt;
617: pseudo->dtctx = ctx;
618: return(0);
619: }
621: /* ----------------------------------------------------------------------------- */
622: /*MC
623: TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping
625: This method solves equations of the form
627: $ F(X,Xdot) = 0
629: for steady state using the iteration
631: $ [G'] S = -F(X,0)
632: $ X += S
634: where
636: $ G(Y) = F(Y,(Y-X)/dt)
638: This is linearly-implicit Euler with the residual always evaluated "at steady
639: state". See note below.
641: Options database keys:
642: + -ts_pseudo_increment <real> - ratio of increase dt
643: - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
645: Level: beginner
647: References:
648: Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
649: C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.
651: Notes:
652: The residual computed by this method includes the transient term (Xdot is computed instead of
653: always being zero), but since the prediction from the last step is always the solution from the
654: last step, on the first Newton iteration we have
656: $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0
658: Therefore, the linear system solved by the first Newton iteration is equivalent to the one
659: described above and in the papers. If the user chooses to perform multiple Newton iterations, the
660: algorithm is no longer the one described in the referenced papers.
662: .seealso: TSCreate(), TS, TSSetType()
664: M*/
667: PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
668: {
669: TS_Pseudo *pseudo;
671: SNES snes;
672: SNESType stype;
675: ts->ops->reset = TSReset_Pseudo;
676: ts->ops->destroy = TSDestroy_Pseudo;
677: ts->ops->view = TSView_Pseudo;
679: ts->ops->setup = TSSetUp_Pseudo;
680: ts->ops->step = TSStep_Pseudo;
681: ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
682: ts->ops->snesfunction = SNESTSFormFunction_Pseudo;
683: ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo;
685: TSGetSNES(ts,&snes);
686: SNESGetType(snes,&stype);
687: if (!stype) {SNESSetType(snes,SNESKSPONLY);}
689: PetscNewLog(ts,&pseudo);
690: ts->data = (void*)pseudo;
692: pseudo->dt_increment = 1.1;
693: pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
694: pseudo->dt = TSPseudoTimeStepDefault;
695: pseudo->fnorm = -1;
697: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);
698: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);
699: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);
700: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);
701: PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);
702: return(0);
703: }
707: /*@C
708: TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
709: Use with TSPseudoSetTimeStep().
711: Collective on TS
713: Input Parameters:
714: . ts - the timestep context
715: . dtctx - unused timestep context
717: Output Parameter:
718: . newdt - the timestep to use for the next step
720: Level: advanced
722: .keywords: timestep, pseudo, default
724: .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
725: @*/
726: PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
727: {
728: TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;
729: PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;
733: VecZeroEntries(pseudo->xdot);
734: TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);
735: VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
736: if (pseudo->fnorm_initial == 0.0) {
737: /* first time through so compute initial function norm */
738: pseudo->fnorm_initial = pseudo->fnorm;
739: fnorm_previous = pseudo->fnorm;
740: }
741: if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step;
742: else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
743: else *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
744: if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
745: pseudo->fnorm_previous = pseudo->fnorm;
746: return(0);
747: }