Actual source code: adaptbasic.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
3: typedef struct {
4: PetscBool always_accept;
5: PetscReal clip[2]; /* admissible decrease/increase factors */
6: PetscReal safety; /* safety factor relative to target error */
7: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
8: Vec Y;
9: } TSAdapt_Basic;
13: static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
14: {
15: TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data;
17: Vec X,Y;
18: PetscReal enorm,hfac_lte,h_lte,safety;
19: PetscInt order,stepno;
22: TSGetTimeStepNumber(ts,&stepno);
23: TSGetSolution(ts,&X);
24: if (!basic->Y) {VecDuplicate(X,&basic->Y);}
25: Y = basic->Y;
26: order = adapt->candidates.order[0];
27: TSEvaluateStep(ts,order-1,Y,NULL);
29: safety = basic->safety;
30: TSErrorNormWRMS(ts,Y,&enorm);
31: if (enorm > 1.) {
32: if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
33: if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
34: PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);
35: *accept = PETSC_TRUE;
36: } else if (basic->always_accept) {
37: PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);
38: *accept = PETSC_TRUE;
39: } else {
40: PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);
41: *accept = PETSC_FALSE;
42: }
43: } else {
44: PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);
45: *accept = PETSC_TRUE;
46: }
48: /* The optimal new step based purely on local truncation error for this step. */
49: hfac_lte = safety * PetscRealPart(PetscPowScalar((PetscScalar)enorm,(PetscReal)(-1./order)));
50: h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
52: *next_sc = 0;
53: *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
54: *wlte = enorm;
55: return(0);
56: }
60: static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
61: {
62: TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data;
66: VecDestroy(&basic->Y);
67: PetscFree(adapt->data);
68: return(0);
69: }
73: static PetscErrorCode TSAdaptSetFromOptions_Basic(TSAdapt adapt)
74: {
75: TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data;
77: PetscInt two;
78: PetscBool set;
81: PetscOptionsHead("Basic adaptive controller options");
82: two = 2;
83: PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);
84: if (set && (two != 2 || basic->clip[0] > basic->clip[1])) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip");
85: PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);
86: PetscOptionsReal("-ts_adapt_basic_reject_safety","Extra safety factor to apply if the last step was rejected","",basic->reject_safety,&basic->reject_safety,NULL);
87: PetscOptionsBool("-ts_adapt_basic_always_accept","Always accept the step regardless of whether local truncation error meets goal","",basic->always_accept,&basic->always_accept,NULL);
88: PetscOptionsTail();
89: return(0);
90: }
94: static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer)
95: {
96: TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data;
98: PetscBool iascii;
101: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
102: if (iascii) {
103: if (basic->always_accept) {PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");}
104: PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);
105: PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);
106: }
107: return(0);
108: }
112: /*MC
113: TSADAPTBASIC - Basic adaptive controller for time stepping
115: Level: intermediate
117: .seealso: TS, TSAdapt, TSSetAdapt()
118: M*/
119: PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
120: {
122: TSAdapt_Basic *a;
125: PetscNewLog(adapt,&a);
126: adapt->data = (void*)a;
127: adapt->ops->choose = TSAdaptChoose_Basic;
128: adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
129: adapt->ops->destroy = TSAdaptDestroy_Basic;
130: adapt->ops->view = TSAdaptView_Basic;
132: a->clip[0] = 0.1;
133: a->clip[1] = 10.;
134: a->safety = 0.9;
135: a->reject_safety = 0.5;
136: a->always_accept = PETSC_FALSE;
137: return(0);
138: }