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