xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
286e171c2SStefano Zampini #include <petscdm.h>
3ae2316d0SEmil Constantinescu 
4ae2316d0SEmil Constantinescu typedef struct {
586e171c2SStefano Zampini   Vec Y;
6726095e4SEmil Constantinescu } TSAdapt_GLEE;
7ae2316d0SEmil Constantinescu 
89371c9d4SSatish Balay static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter) {
9657c1e31SEmil Constantinescu   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;
108a175baeSEmil Constantinescu   Vec           X, Y, E;
11df3bac28SEmil Constantinescu   PetscReal     enorm, enorma, enormr, hfac_lte, hfac_ltea, hfac_lter, h_lte, safety;
1280275a0aSLisandro Dalcin   PetscInt      order;
1386e171c2SStefano Zampini   PetscBool     bGTEMethod;
14ae2316d0SEmil Constantinescu 
15ae2316d0SEmil Constantinescu   PetscFunctionBegin;
165e4ed32fSEmil Constantinescu   *next_sc = 0; /* Reuse the same order scheme */
171917a363SLisandro Dalcin   safety   = adapt->safety;
189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSGLEE, &bGTEMethod));
19df3bac28SEmil Constantinescu   order = adapt->candidates.order[0];
20df3bac28SEmil Constantinescu 
21df3bac28SEmil Constantinescu   if (bGTEMethod) { /* the method is of GLEE type */
2286e171c2SStefano Zampini     DM dm;
2386e171c2SStefano Zampini 
249566063dSJacob Faibussowitsch     PetscCall(TSGetSolution(ts, &X));
251c167fc2SEmil Constantinescu     if (!glee->Y && adapt->glee_use_local) {
269566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X, &glee->Y)); /*create vector to store previous step global error*/
279566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(glee->Y));   /*set error to zero on the first step - may not work if error is not zero initially*/
281c167fc2SEmil Constantinescu     }
299566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
309566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &E));
319566063dSJacob Faibussowitsch     PetscCall(TSGetTimeError(ts, 0, &E));
321c167fc2SEmil Constantinescu 
339566063dSJacob Faibussowitsch     if (adapt->glee_use_local) PetscCall(VecAXPY(E, -1.0, glee->Y)); /* local error = current error - previous step error */
341c167fc2SEmil Constantinescu 
351c167fc2SEmil Constantinescu     /* this should be called with the solution at the beginning of the step too*/
369566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedENorm(ts, E, X, X, adapt->wnormtype, &enorm, &enorma, &enormr));
379566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &E));
380a01e1b2SEmil Constantinescu   } else {
39df3bac28SEmil Constantinescu     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
409566063dSJacob Faibussowitsch     PetscCall(TSGetSolution(ts, &X));
419566063dSJacob Faibussowitsch     if (!glee->Y) PetscCall(VecDuplicate(X, &glee->Y));
42657c1e31SEmil Constantinescu     Y = glee->Y;
439566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL));
449566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, adapt->wnormtype, &enorm, &enorma, &enormr));
450a01e1b2SEmil Constantinescu   }
465e4ed32fSEmil Constantinescu 
475e4ed32fSEmil Constantinescu   if (enorm < 0) {
485e4ed32fSEmil Constantinescu     *accept = PETSC_TRUE;
495e4ed32fSEmil Constantinescu     *next_h = h;  /* Reuse the old step */
50df3bac28SEmil Constantinescu     *wlte   = -1; /* Weighted error was not evaluated */
51df3bac28SEmil Constantinescu     *wltea  = -1; /* Weighted absolute error was not evaluated */
52df3bac28SEmil Constantinescu     *wlter  = -1; /* Weighted relative error was not evaluated */
535e4ed32fSEmil Constantinescu     PetscFunctionReturn(0);
545e4ed32fSEmil Constantinescu   }
555e4ed32fSEmil Constantinescu 
56df3bac28SEmil Constantinescu   if (enorm > 1. || enorma > 1. || enormr > 1.) {
571917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
58ae2316d0SEmil Constantinescu     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON) * adapt->dt_min) {
599566063dSJacob Faibussowitsch       PetscCall(PetscInfo(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));
60ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
61bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(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));
63ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
64ae2316d0SEmil Constantinescu     } else {
659566063dSJacob Faibussowitsch       PetscCall(PetscInfo(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));
66ae2316d0SEmil Constantinescu       *accept = PETSC_FALSE;
67ae2316d0SEmil Constantinescu     }
68ae2316d0SEmil Constantinescu   } else {
699566063dSJacob Faibussowitsch     PetscCall(PetscInfo(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));
70ae2316d0SEmil Constantinescu     *accept = PETSC_TRUE;
71ae2316d0SEmil Constantinescu   }
72ae2316d0SEmil Constantinescu 
73df3bac28SEmil Constantinescu   if (bGTEMethod) {
741c167fc2SEmil Constantinescu     if (*accept == PETSC_TRUE && adapt->glee_use_local) {
751c167fc2SEmil Constantinescu       /* If step is accepted, then overwrite previous step error with the current error to be used on the next step */
761c167fc2SEmil Constantinescu       /* WARNING: if the adapters are composable, then the accept test will not be reliable*/
779566063dSJacob Faibussowitsch       PetscCall(TSGetTimeError(ts, 0, &glee->Y));
781c167fc2SEmil Constantinescu     }
791c167fc2SEmil Constantinescu 
80df3bac28SEmil Constantinescu     /* The optimal new step based on the current global truncation error. */
81df3bac28SEmil Constantinescu     if (enorm > 0) {
82df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
831c167fc2SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(1. / enorma, ((PetscReal)1) / (order + 1));
84df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
851c167fc2SEmil Constantinescu       hfac_lter = safety * PetscPowReal(1. / enormr, ((PetscReal)1) / (order + 1));
86df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
87df3bac28SEmil Constantinescu       hfac_lte  = PetscMin(hfac_ltea, hfac_lter);
88df3bac28SEmil Constantinescu     } else {
89ae2316d0SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
90df3bac28SEmil Constantinescu     }
911917a363SLisandro Dalcin     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
92ae2316d0SEmil Constantinescu     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
93df3bac28SEmil Constantinescu   } else {
94df3bac28SEmil Constantinescu     /* The optimal new step based purely on local truncation error for this step. */
95df3bac28SEmil Constantinescu     if (enorm > 0) {
96df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
97df3bac28SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(enorma, ((PetscReal)-1) / order);
98df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
99df3bac28SEmil Constantinescu       hfac_lter = safety * PetscPowReal(enormr, ((PetscReal)-1) / order);
100df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
101df3bac28SEmil Constantinescu       hfac_lte  = PetscMin(hfac_ltea, hfac_lter);
102df3bac28SEmil Constantinescu     } else {
103df3bac28SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
104df3bac28SEmil Constantinescu     }
1051917a363SLisandro Dalcin     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
106df3bac28SEmil Constantinescu     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
107df3bac28SEmil Constantinescu   }
108ae2316d0SEmil Constantinescu   *wlte  = enorm;
1095e4ed32fSEmil Constantinescu   *wltea = enorma;
1105e4ed32fSEmil Constantinescu   *wlter = enormr;
111ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
112ae2316d0SEmil Constantinescu }
113ae2316d0SEmil Constantinescu 
1149371c9d4SSatish Balay static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) {
115657c1e31SEmil Constantinescu   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;
116ae2316d0SEmil Constantinescu 
117ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->Y));
119ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
120ae2316d0SEmil Constantinescu }
121ae2316d0SEmil Constantinescu 
1229371c9d4SSatish Balay static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) {
123ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(TSAdaptReset_GLEE(adapt));
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree(adapt->data));
126ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
127ae2316d0SEmil Constantinescu }
128ae2316d0SEmil Constantinescu 
129ae2316d0SEmil Constantinescu /*MC
130657c1e31SEmil Constantinescu    TSADAPTGLEE - GLEE adaptive controller for time stepping
131ae2316d0SEmil Constantinescu 
132ae2316d0SEmil Constantinescu    Level: intermediate
133ae2316d0SEmil Constantinescu 
134db781477SPatrick Sanan .seealso: `TS`, `TSAdapt`, `TSGetAdapt()`
135ae2316d0SEmil Constantinescu M*/
1369371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) {
1371917a363SLisandro Dalcin   TSAdapt_GLEE *glee;
138ae2316d0SEmil Constantinescu 
139ae2316d0SEmil Constantinescu   PetscFunctionBegin;
140*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&glee));
1411917a363SLisandro Dalcin   adapt->data         = (void *)glee;
142726095e4SEmil Constantinescu   adapt->ops->choose  = TSAdaptChoose_GLEE;
1431917a363SLisandro Dalcin   adapt->ops->reset   = TSAdaptReset_GLEE;
144726095e4SEmil Constantinescu   adapt->ops->destroy = TSAdaptDestroy_GLEE;
145ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
146ae2316d0SEmil Constantinescu }
147