xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision f2300f31d01ffc99ca7a05a43941a14c1c1d9061)
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 
TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt * next_sc,PetscReal * next_h,PetscBool * accept,PetscReal * wlte,PetscReal * wltea,PetscReal * wlter)8d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
9d71ae5a4SJacob Faibussowitsch {
10657c1e31SEmil Constantinescu   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;
118a175baeSEmil Constantinescu   Vec           X, Y, E;
12df3bac28SEmil Constantinescu   PetscReal     enorm, enorma, enormr, hfac_lte, hfac_ltea, hfac_lter, h_lte, safety;
1380275a0aSLisandro Dalcin   PetscInt      order;
1486e171c2SStefano Zampini   PetscBool     bGTEMethod;
15ae2316d0SEmil Constantinescu 
16ae2316d0SEmil Constantinescu   PetscFunctionBegin;
175e4ed32fSEmil Constantinescu   *next_sc = 0; /* Reuse the same order scheme */
181917a363SLisandro Dalcin   safety   = adapt->safety;
199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSGLEE, &bGTEMethod));
20df3bac28SEmil Constantinescu   order = adapt->candidates.order[0];
21df3bac28SEmil Constantinescu 
22df3bac28SEmil Constantinescu   if (bGTEMethod) { /* the method is of GLEE type */
2386e171c2SStefano Zampini     DM dm;
2486e171c2SStefano Zampini 
259566063dSJacob Faibussowitsch     PetscCall(TSGetSolution(ts, &X));
261c167fc2SEmil Constantinescu     if (!glee->Y && adapt->glee_use_local) {
279566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X, &glee->Y)); /*create vector to store previous step global error*/
289566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(glee->Y));   /*set error to zero on the first step - may not work if error is not zero initially*/
291c167fc2SEmil Constantinescu     }
309566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
319566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &E));
329566063dSJacob Faibussowitsch     PetscCall(TSGetTimeError(ts, 0, &E));
331c167fc2SEmil Constantinescu 
349566063dSJacob Faibussowitsch     if (adapt->glee_use_local) PetscCall(VecAXPY(E, -1.0, glee->Y)); /* local error = current error - previous step error */
351c167fc2SEmil Constantinescu 
361c167fc2SEmil Constantinescu     /* this should be called with the solution at the beginning of the step too*/
379566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedENorm(ts, E, X, X, adapt->wnormtype, &enorm, &enorma, &enormr));
389566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &E));
390a01e1b2SEmil Constantinescu   } else {
40df3bac28SEmil Constantinescu     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
419566063dSJacob Faibussowitsch     PetscCall(TSGetSolution(ts, &X));
429566063dSJacob Faibussowitsch     if (!glee->Y) PetscCall(VecDuplicate(X, &glee->Y));
43657c1e31SEmil Constantinescu     Y = glee->Y;
449566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL));
459566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, adapt->wnormtype, &enorm, &enorma, &enormr));
460a01e1b2SEmil Constantinescu   }
475e4ed32fSEmil Constantinescu 
485e4ed32fSEmil Constantinescu   if (enorm < 0) {
495e4ed32fSEmil Constantinescu     *accept = PETSC_TRUE;
505e4ed32fSEmil Constantinescu     *next_h = h;  /* Reuse the old step */
51df3bac28SEmil Constantinescu     *wlte   = -1; /* Weighted error was not evaluated */
52df3bac28SEmil Constantinescu     *wltea  = -1; /* Weighted absolute error was not evaluated */
53df3bac28SEmil Constantinescu     *wlter  = -1; /* Weighted relative error was not evaluated */
543ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
555e4ed32fSEmil Constantinescu   }
565e4ed32fSEmil Constantinescu 
57df3bac28SEmil Constantinescu   if (enorm > 1. || enorma > 1. || enormr > 1.) {
581917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
59ae2316d0SEmil Constantinescu     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON) * adapt->dt_min) {
609566063dSJacob 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));
61ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
62bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
639566063dSJacob 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));
64ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
65ae2316d0SEmil Constantinescu     } else {
669566063dSJacob 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));
67ae2316d0SEmil Constantinescu       *accept = PETSC_FALSE;
68ae2316d0SEmil Constantinescu     }
69ae2316d0SEmil Constantinescu   } else {
709566063dSJacob 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));
71ae2316d0SEmil Constantinescu     *accept = PETSC_TRUE;
72ae2316d0SEmil Constantinescu   }
73ae2316d0SEmil Constantinescu 
74df3bac28SEmil Constantinescu   if (bGTEMethod) {
751c167fc2SEmil Constantinescu     if (*accept == PETSC_TRUE && adapt->glee_use_local) {
761c167fc2SEmil Constantinescu       /* If step is accepted, then overwrite previous step error with the current error to be used on the next step */
771c167fc2SEmil Constantinescu       /* WARNING: if the adapters are composable, then the accept test will not be reliable*/
789566063dSJacob Faibussowitsch       PetscCall(TSGetTimeError(ts, 0, &glee->Y));
791c167fc2SEmil Constantinescu     }
801c167fc2SEmil Constantinescu 
81df3bac28SEmil Constantinescu     /* The optimal new step based on the current global truncation error. */
82df3bac28SEmil Constantinescu     if (enorm > 0) {
83df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
841c167fc2SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(1. / enorma, ((PetscReal)1) / (order + 1));
85df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
861c167fc2SEmil Constantinescu       hfac_lter = safety * PetscPowReal(1. / enormr, ((PetscReal)1) / (order + 1));
87df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
88df3bac28SEmil Constantinescu       hfac_lte = PetscMin(hfac_ltea, hfac_lter);
89df3bac28SEmil Constantinescu     } else {
90ae2316d0SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
91df3bac28SEmil Constantinescu     }
921917a363SLisandro Dalcin     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
93ae2316d0SEmil Constantinescu     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
94df3bac28SEmil Constantinescu   } else {
95df3bac28SEmil Constantinescu     /* The optimal new step based purely on local truncation error for this step. */
96df3bac28SEmil Constantinescu     if (enorm > 0) {
97df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
98df3bac28SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(enorma, ((PetscReal)-1) / order);
99df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
100df3bac28SEmil Constantinescu       hfac_lter = safety * PetscPowReal(enormr, ((PetscReal)-1) / order);
101df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
102df3bac28SEmil Constantinescu       hfac_lte = PetscMin(hfac_ltea, hfac_lter);
103df3bac28SEmil Constantinescu     } else {
104df3bac28SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
105df3bac28SEmil Constantinescu     }
1061917a363SLisandro Dalcin     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
107df3bac28SEmil Constantinescu     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
108df3bac28SEmil Constantinescu   }
109ae2316d0SEmil Constantinescu   *wlte  = enorm;
1105e4ed32fSEmil Constantinescu   *wltea = enorma;
1115e4ed32fSEmil Constantinescu   *wlter = enormr;
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113ae2316d0SEmil Constantinescu }
114ae2316d0SEmil Constantinescu 
TSAdaptReset_GLEE(TSAdapt adapt)115d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
116d71ae5a4SJacob Faibussowitsch {
117657c1e31SEmil Constantinescu   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;
118ae2316d0SEmil Constantinescu 
119ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->Y));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122ae2316d0SEmil Constantinescu }
123ae2316d0SEmil Constantinescu 
TSAdaptDestroy_GLEE(TSAdapt adapt)124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt)
125d71ae5a4SJacob Faibussowitsch {
126ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1279566063dSJacob Faibussowitsch   PetscCall(TSAdaptReset_GLEE(adapt));
1289566063dSJacob Faibussowitsch   PetscCall(PetscFree(adapt->data));
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130ae2316d0SEmil Constantinescu }
131ae2316d0SEmil Constantinescu 
132ae2316d0SEmil Constantinescu /*MC
133657c1e31SEmil Constantinescu    TSADAPTGLEE - GLEE adaptive controller for time stepping
134ae2316d0SEmil Constantinescu 
135ae2316d0SEmil Constantinescu    Level: intermediate
136ae2316d0SEmil Constantinescu 
137*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`
138ae2316d0SEmil Constantinescu M*/
TSAdaptCreate_GLEE(TSAdapt adapt)139d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt)
140d71ae5a4SJacob Faibussowitsch {
1411917a363SLisandro Dalcin   TSAdapt_GLEE *glee;
142ae2316d0SEmil Constantinescu 
143ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1444dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&glee));
1451917a363SLisandro Dalcin   adapt->data         = (void *)glee;
146726095e4SEmil Constantinescu   adapt->ops->choose  = TSAdaptChoose_GLEE;
1471917a363SLisandro Dalcin   adapt->ops->reset   = TSAdaptReset_GLEE;
148726095e4SEmil Constantinescu   adapt->ops->destroy = TSAdaptDestroy_GLEE;
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150ae2316d0SEmil Constantinescu }
151