xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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       X;
9   Vec       Y;
10   Vec       E;
11 } TSAdapt_GLEE;
12 
13 static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
14 {
15   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
16   TSType         time_scheme;      /* Type of time-integration scheme        */
17   PetscErrorCode ierr;
18   Vec            X,Y,E;
19   PetscReal      enorm,enorma,enormr,hfac_lte,hfac_ltea,hfac_lter,h_lte,safety;
20   PetscInt       order,stepno;
21   PetscBool      bGTEMethod=PETSC_FALSE;
22 
23   PetscFunctionBegin;
24 
25   *next_sc = 0; /* Reuse the same order scheme */
26   safety = glee->safety;
27   ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr);
28   ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr);
29   if (!strcmp(time_scheme,TSGLEE)) bGTEMethod=PETSC_TRUE;
30   order = adapt->candidates.order[0];
31 
32   if (bGTEMethod){/* the method is of GLEE type */
33     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
34     if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);}
35     E     = glee->E;
36     ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr);
37     /* this should be called with Y (the solution at the beginning of the step)*/
38     ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
39   } else {
40     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
41     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
42     if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);}
43     Y     = glee->Y;
44     ierr  = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr);
45     ierr  = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
46   }
47 
48   if (enorm < 0) {
49     *accept  = PETSC_TRUE;
50     *next_h  = h;            /* Reuse the old step */
51     *wlte    = -1;           /* Weighted error was not evaluated */
52     *wltea   = -1;           /* Weighted absolute error was not evaluated */
53     *wlter   = -1;           /* Weighted relative error was not evaluated */
54     PetscFunctionReturn(0);
55   }
56 
57   if (enorm > 1. || enorma > 1. || enormr > 1.) {
58     if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */
59     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
60       ierr    = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting because step size %g is at minimum\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr);
61       *accept = PETSC_TRUE;
62     } else if (glee->always_accept) {
63       ierr    = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting step of size %g because always_accept is set\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr);
64       *accept = PETSC_TRUE;
65     } else {
66       ierr    = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], rejecting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr);
67       *accept = PETSC_FALSE;
68     }
69   } else {
70     ierr    = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative] [%g, %g, %g], accepting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr);
71     *accept = PETSC_TRUE;
72   }
73 
74   if (bGTEMethod){
75     /* The optimal new step based on the current global truncation error. */
76     if (enorm > 0) {
77       /* factor based on the absolute tolerance */
78       hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/order);
79       /* factor based on the relative tolerance */
80       hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/order);
81       /* pick the minimum time step among the relative and absolute tolerances */
82       hfac_lte  = PetscMin(hfac_ltea,hfac_lter);
83     } else {
84       hfac_lte = safety * PETSC_INFINITY;
85     }
86     h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]);
87     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
88   } else {
89     /* The optimal new step based purely on local truncation error for this step. */
90     if (enorm > 0) {
91       /* factor based on the absolute tolerance */
92       hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order);
93       /* factor based on the relative tolerance */
94       hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order);
95       /* pick the minimum time step among the relative and absolute tolerances */
96       hfac_lte  = PetscMin(hfac_ltea,hfac_lter);
97     } else {
98       hfac_lte = safety * PETSC_INFINITY;
99     }
100     h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]);
101     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
102   }
103   *wlte   = enorm;
104   *wltea  = enorma;
105   *wlter  = enormr;
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
110 {
111   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   ierr = VecDestroy(&glee->Y);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt)
120 {
121   PetscErrorCode ierr;
122 
123   PetscFunctionBegin;
124   ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr);
125   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
130 {
131   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
132   PetscErrorCode ierr;
133   PetscInt       two;
134   PetscBool      set;
135 
136   PetscFunctionBegin;
137   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr);
138   two  = 2;
139   ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr);
140   if (set && (two != 2 || glee->clip[0] > glee->clip[1])) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_glee_clip");
141   ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr);
142   ierr = PetscOptionsReal("-ts_adapt_glee_reject_safety","Extra safety factor to apply if the last step was rejected","",glee->reject_safety,&glee->reject_safety,NULL);CHKERRQ(ierr);
143   ierr = PetscOptionsBool("-ts_adapt_glee_always_accept","Always accept the step regardless of whether local truncation error meets goal","",glee->always_accept,&glee->always_accept,NULL);CHKERRQ(ierr);
144   ierr = PetscOptionsTail();CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer)
149 {
150   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
151   PetscErrorCode ierr;
152   PetscBool      iascii;
153 
154   PetscFunctionBegin;
155   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
156   if (iascii) {
157     if (glee->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: always accepting steps\n");CHKERRQ(ierr);}
158     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr);
159     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr);
160   }
161   PetscFunctionReturn(0);
162 }
163 
164 /*MC
165    TSADAPTGLEE - GLEE adaptive controller for time stepping
166 
167    Level: intermediate
168 
169 .seealso: TS, TSAdapt, TSSetAdapt()
170 M*/
171 PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt)
172 {
173   PetscErrorCode ierr;
174   TSAdapt_GLEE  *a;
175 
176   PetscFunctionBegin;
177   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
178   adapt->data                = (void*)a;
179   adapt->ops->choose         = TSAdaptChoose_GLEE;
180   adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE;
181   adapt->ops->destroy        = TSAdaptDestroy_GLEE;
182   adapt->ops->view           = TSAdaptView_GLEE;
183 
184   a->clip[0]       = 0.1;
185   a->clip[1]       = 10.;
186   a->safety        = 0.9;
187   a->reject_safety = 0.5;
188   a->always_accept = PETSC_FALSE;
189   PetscFunctionReturn(0);
190 }
191