xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision 966be33a19c9230d4aa438248a644248d45cc287)
1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2 
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;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "TSAdaptChoose_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;
16   PetscInt       order  = PETSC_DECIDE;
17   PetscReal      enorm  = -1;
18   PetscReal      safety = basic->safety;
19   PetscReal      hfac_lte,h_lte;
20   PetscErrorCode ierr;
21 
22   PetscFunctionBegin;
23   *next_sc = 0; /* Reuse the same order scheme */
24 
25   if (ts->ops->evaluatewlte) {
26     ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr);
27     if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order);
28   } else if (ts->ops->evaluatestep) {
29     if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered");
30     if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
31     if (!basic->Y) {ierr = VecDuplicate(ts->vec_sol,&basic->Y);CHKERRQ(ierr);}
32     order = adapt->candidates.order[0];
33     ierr = TSEvaluateStep(ts,order-1,basic->Y,NULL);CHKERRQ(ierr);
34     ierr = TSErrorWeightedNorm(ts,ts->vec_sol,basic->Y,adapt->wnormtype,&enorm);CHKERRQ(ierr);
35   }
36 
37   if (enorm < 0) {
38     *accept  = PETSC_TRUE;
39     *next_h  = h;            /* Reuse the old step */
40     *wlte    = -1;           /* Weighted local truncation error was not evaluated */
41     PetscFunctionReturn(0);
42   }
43 
44   /* Determine whether the step is accepted of rejected */
45   if (enorm > 1) {
46     if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
47     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
48       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr);
49       *accept = PETSC_TRUE;
50     } else if (basic->always_accept) {
51       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);CHKERRQ(ierr);
52       *accept = PETSC_TRUE;
53     } else {
54       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
55       *accept = PETSC_FALSE;
56     }
57   } else {
58     ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
59     *accept = PETSC_TRUE;
60   }
61 
62   /* The optimal new step based purely on local truncation error for this step. */
63   if (enorm > 0)
64     hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order);
65   else
66     hfac_lte = safety * PETSC_INFINITY;
67   h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
68 
69   *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
70   *wlte   = enorm;
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "TSAdaptReset_Basic"
76 static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt)
77 {
78   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "TSAdaptDestroy_Basic"
88 static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
89 {
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr);
94   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 
98 #undef __FUNCT__
99 #define __FUNCT__ "TSAdaptSetFromOptions_Basic"
100 static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
101 {
102   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
103   PetscErrorCode ierr;
104   PetscInt       two;
105   PetscBool      set;
106 
107   PetscFunctionBegin;
108   ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr);
109   two  = 2;
110   ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase factor in step size","TSAdaptBasicSetClip",basic->clip,&two,&set);CHKERRQ(ierr);
111   if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip");
112   if (set) {ierr = TSAdaptBasicSetClip(adapt,basic->clip[0],basic->clip[1]);CHKERRQ(ierr);}
113   ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr);
114   ierr = PetscOptionsReal("-ts_adapt_basic_reject_safety","Extra safety factor to apply if the last step was rejected","",basic->reject_safety,&basic->reject_safety,NULL);CHKERRQ(ierr);
115   ierr = 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);CHKERRQ(ierr);
116   ierr = PetscOptionsTail();CHKERRQ(ierr);
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "TSAdaptView_Basic"
122 static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer)
123 {
124   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
125   PetscErrorCode ierr;
126   PetscBool      iascii;
127 
128   PetscFunctionBegin;
129   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
130   if (iascii) {
131     if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  Basic: always accepting steps\n");CHKERRQ(ierr);}
132     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr);
133     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr);
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "TSAdaptCreate_Basic"
140 /*MC
141    TSADAPTBASIC - Basic adaptive controller for time stepping
142 
143    Level: intermediate
144 
145 .seealso: TS, TSAdapt, TSSetAdapt()
146 M*/
147 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
148 {
149   PetscErrorCode ierr;
150   TSAdapt_Basic  *a;
151 
152   PetscFunctionBegin;
153   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
154   adapt->data                = (void*)a;
155   adapt->ops->choose         = TSAdaptChoose_Basic;
156   adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
157   adapt->ops->destroy        = TSAdaptDestroy_Basic;
158   adapt->ops->view           = TSAdaptView_Basic;
159 
160   a->clip[0]       = 0.1;
161   a->clip[1]       = 10.;
162   a->safety        = 0.9;
163   a->reject_safety = 0.5;
164   a->always_accept = PETSC_FALSE;
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNCT__
169 #define __FUNCT__ "TSAdaptBasicSetClip"
170 /*@
171    TSAdaptBasicSetClip - Sets the admissible decrease/increase factor in step size
172 
173    Collective on TSAdapt
174 
175    Input Arguments:
176 +  adapt - adaptive controller context
177 .  low - admissible decrease factor
178 +  high - admissible increase factor
179 
180    Level: intermediate
181 
182 .seealso: TSAdaptChoose()
183 @*/
184 PetscErrorCode TSAdaptBasicSetClip(TSAdapt adapt,PetscReal low,PetscReal high)
185 {
186   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
187   PetscBool      isbasic;
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
192   if (low  != PETSC_DEFAULT && low  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low);
193   if (low  != PETSC_DEFAULT && low  > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low);
194   if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high);
195   ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr);
196   if (isbasic) {
197     if (low  != PETSC_DEFAULT) basic->clip[0] = low;
198     if (high != PETSC_DEFAULT) basic->clip[1] = high;
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "TSAdaptBasicGetClip"
205 /*@
206    TSAdaptBasicGetClip - Gets the admissible decrease/increase factor in step size
207 
208    Collective on TSAdapt
209 
210    Input Arguments:
211 .  adapt - adaptive controller context
212 
213    Ouput Arguments:
214 +  low - optional, admissible decrease factor
215 -  high - optional, admissible increase factor
216 
217    Level: intermediate
218 
219 .seealso: TSAdaptChoose()
220 @*/
221 PetscErrorCode TSAdaptBasicGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high)
222 {
223   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
224   PetscBool      isbasic;
225   PetscErrorCode ierr;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
229   if (low)  PetscValidRealPointer(low,2);
230   if (high) PetscValidRealPointer(high,3);
231   ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr);
232   if (isbasic) {
233     if (low)  *low  = basic->clip[0];
234     if (high) *high = basic->clip[1];
235   } else {
236     if (low)  *low  = 0;
237     if (high) *high = PETSC_MAX_REAL;
238   }
239   PetscFunctionReturn(0);
240 }
241