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