xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision 0a01e1b289b39f06e756d468ae9070cb2b3504ae)
1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2ae2316d0SEmil Constantinescu 
3ae2316d0SEmil Constantinescu typedef struct {
4ae2316d0SEmil Constantinescu   PetscBool always_accept;
5ae2316d0SEmil Constantinescu   PetscReal clip[2];            /* admissible decrease/increase factors */
6ae2316d0SEmil Constantinescu   PetscReal safety;             /* safety factor relative to target error */
7ae2316d0SEmil Constantinescu   PetscReal reject_safety;      /* extra safety factor if the last step was rejected */
8ae2316d0SEmil Constantinescu   Vec       Y;
9726095e4SEmil Constantinescu } TSAdapt_GLEE;
10ae2316d0SEmil Constantinescu 
11ae2316d0SEmil Constantinescu #undef __FUNCT__
12726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptChoose_GLEE"
13726095e4SEmil Constantinescu static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
14ae2316d0SEmil Constantinescu {
15726095e4SEmil Constantinescu   TSAdapt_GLEE  *basic = (TSAdapt_GLEE*)adapt->data;
16*0a01e1b2SEmil Constantinescu   TSType         time_scheme;      /* Type of time-integration scheme        */
17ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
18ae2316d0SEmil Constantinescu   Vec            X,Y;
197453f775SEmil Constantinescu   PetscReal      enorm,enorma,enormr,hfac_lte,h_lte,safety;
20ae2316d0SEmil Constantinescu   PetscInt       order,stepno;
21ae2316d0SEmil Constantinescu 
22ae2316d0SEmil Constantinescu   PetscFunctionBegin;
23ae2316d0SEmil Constantinescu   ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr);
24*0a01e1b2SEmil Constantinescu 
25*0a01e1b2SEmil Constantinescu   safety = basic->safety;
26*0a01e1b2SEmil Constantinescu   ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr);
27*0a01e1b2SEmil Constantinescu   if (!strcmp(time_scheme,TSGLEE)){
28*0a01e1b2SEmil Constantinescu     /* the method is of GLEE type */
29*0a01e1b2SEmil Constantinescu     ierr = TSGetTimeError(ts,-1,&X);CHKERRQ(ierr);
30*0a01e1b2SEmil Constantinescu     ierr = TSGetTimeError(ts, 0,&Y);CHKERRQ(ierr);
31*0a01e1b2SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
32*0a01e1b2SEmil Constantinescu   } else {
33*0a01e1b2SEmil Constantinescu     /* the method is NOT of GLEE type */
34ae2316d0SEmil Constantinescu     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
35ae2316d0SEmil Constantinescu     if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);}
36ae2316d0SEmil Constantinescu     Y     = basic->Y;
37ae2316d0SEmil Constantinescu     order = adapt->candidates.order[0];
38ae2316d0SEmil Constantinescu     ierr  = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr);
397453f775SEmil Constantinescu     ierr  = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
40*0a01e1b2SEmil Constantinescu   }
41ae2316d0SEmil Constantinescu   if (enorm > 1.) {
42ae2316d0SEmil Constantinescu     if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
43ae2316d0SEmil Constantinescu     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
44ae2316d0SEmil Constantinescu       ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr);
45ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
46ae2316d0SEmil Constantinescu     } else if (basic->always_accept) {
47ae2316d0SEmil Constantinescu       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);
48ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
49ae2316d0SEmil Constantinescu     } else {
50ae2316d0SEmil Constantinescu       ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
51ae2316d0SEmil Constantinescu       *accept = PETSC_FALSE;
52ae2316d0SEmil Constantinescu     }
53ae2316d0SEmil Constantinescu   } else {
54ae2316d0SEmil Constantinescu     ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
55ae2316d0SEmil Constantinescu     *accept = PETSC_TRUE;
56ae2316d0SEmil Constantinescu   }
57ae2316d0SEmil Constantinescu 
58ae2316d0SEmil Constantinescu   /* The optimal new step based purely on local truncation error for this step. */
59ae2316d0SEmil Constantinescu   if (enorm == 0.0) {
60ae2316d0SEmil Constantinescu     hfac_lte = safety * PETSC_INFINITY;
61ae2316d0SEmil Constantinescu   } else {
62ae2316d0SEmil Constantinescu     hfac_lte = safety * PetscPowReal(enorm,-1./order);
63ae2316d0SEmil Constantinescu   }
64ae2316d0SEmil Constantinescu   h_lte    = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
65ae2316d0SEmil Constantinescu 
66ae2316d0SEmil Constantinescu   *next_sc = 0;
67ae2316d0SEmil Constantinescu   *next_h  = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
68ae2316d0SEmil Constantinescu   *wlte    = enorm;
69ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
70ae2316d0SEmil Constantinescu }
71ae2316d0SEmil Constantinescu 
72ae2316d0SEmil Constantinescu #undef __FUNCT__
73726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptReset_GLEE"
74726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
75ae2316d0SEmil Constantinescu {
76726095e4SEmil Constantinescu   TSAdapt_GLEE  *basic = (TSAdapt_GLEE*)adapt->data;
77ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
78ae2316d0SEmil Constantinescu 
79ae2316d0SEmil Constantinescu   PetscFunctionBegin;
80ae2316d0SEmil Constantinescu   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
81ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
82ae2316d0SEmil Constantinescu }
83ae2316d0SEmil Constantinescu 
84ae2316d0SEmil Constantinescu #undef __FUNCT__
85726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptDestroy_GLEE"
86726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt)
87ae2316d0SEmil Constantinescu {
88ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
89ae2316d0SEmil Constantinescu 
90ae2316d0SEmil Constantinescu   PetscFunctionBegin;
91726095e4SEmil Constantinescu   ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr);
92ae2316d0SEmil Constantinescu   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
93ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
94ae2316d0SEmil Constantinescu }
95ae2316d0SEmil Constantinescu 
96ae2316d0SEmil Constantinescu #undef __FUNCT__
97726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptSetFromOptions_GLEE"
98726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
99ae2316d0SEmil Constantinescu {
100726095e4SEmil Constantinescu   TSAdapt_GLEE  *basic = (TSAdapt_GLEE*)adapt->data;
101ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
102ae2316d0SEmil Constantinescu   PetscInt       two;
103ae2316d0SEmil Constantinescu   PetscBool      set;
104ae2316d0SEmil Constantinescu 
105ae2316d0SEmil Constantinescu   PetscFunctionBegin;
106726095e4SEmil Constantinescu   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr);
107ae2316d0SEmil Constantinescu   two  = 2;
108ae2316d0SEmil Constantinescu   ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr);
109ae2316d0SEmil Constantinescu   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");
110ae2316d0SEmil Constantinescu   ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr);
111ae2316d0SEmil Constantinescu   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);
112ae2316d0SEmil Constantinescu   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);
113ae2316d0SEmil Constantinescu   ierr = PetscOptionsTail();CHKERRQ(ierr);
114ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
115ae2316d0SEmil Constantinescu }
116ae2316d0SEmil Constantinescu 
117ae2316d0SEmil Constantinescu #undef __FUNCT__
118726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptView_GLEE"
119726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer)
120ae2316d0SEmil Constantinescu {
121726095e4SEmil Constantinescu   TSAdapt_GLEE  *basic = (TSAdapt_GLEE*)adapt->data;
122ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
123ae2316d0SEmil Constantinescu   PetscBool      iascii;
124ae2316d0SEmil Constantinescu 
125ae2316d0SEmil Constantinescu   PetscFunctionBegin;
126ae2316d0SEmil Constantinescu   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
127ae2316d0SEmil Constantinescu   if (iascii) {
128726095e4SEmil Constantinescu     if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: always accepting steps\n");CHKERRQ(ierr);}
129726095e4SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr);
130726095e4SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr);
131ae2316d0SEmil Constantinescu   }
132ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
133ae2316d0SEmil Constantinescu }
134ae2316d0SEmil Constantinescu 
135ae2316d0SEmil Constantinescu #undef __FUNCT__
136726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptCreate_GLEE"
137ae2316d0SEmil Constantinescu /*MC
138726095e4SEmil Constantinescu    TSADAPTBASIC - GLEE adaptive controller for time stepping
139ae2316d0SEmil Constantinescu 
140ae2316d0SEmil Constantinescu    Level: intermediate
141ae2316d0SEmil Constantinescu 
142ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt()
143ae2316d0SEmil Constantinescu M*/
144726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt)
145ae2316d0SEmil Constantinescu {
146ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
147726095e4SEmil Constantinescu   TSAdapt_GLEE  *a;
148ae2316d0SEmil Constantinescu 
149ae2316d0SEmil Constantinescu   PetscFunctionBegin;
150ae2316d0SEmil Constantinescu   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
151ae2316d0SEmil Constantinescu   adapt->data                = (void*)a;
152726095e4SEmil Constantinescu   adapt->ops->choose         = TSAdaptChoose_GLEE;
153726095e4SEmil Constantinescu   adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE;
154726095e4SEmil Constantinescu   adapt->ops->destroy        = TSAdaptDestroy_GLEE;
155726095e4SEmil Constantinescu   adapt->ops->view           = TSAdaptView_GLEE;
156ae2316d0SEmil Constantinescu 
157ae2316d0SEmil Constantinescu   a->clip[0]       = 0.1;
158ae2316d0SEmil Constantinescu   a->clip[1]       = 10.;
159ae2316d0SEmil Constantinescu   a->safety        = 0.9;
160ae2316d0SEmil Constantinescu   a->reject_safety = 0.5;
161ae2316d0SEmil Constantinescu   a->always_accept = PETSC_FALSE;
162ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
163ae2316d0SEmil Constantinescu }
164