xref: /petsc/src/ts/impls/implicit/glle/glleadapt.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
13d177a5cSEmil Constantinescu 
23d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I  "petscts.h" I*/
33d177a5cSEmil Constantinescu 
43d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAdaptList;
53d177a5cSEmil Constantinescu static PetscBool         TSGLLEAdaptPackageInitialized;
63d177a5cSEmil Constantinescu static PetscBool         TSGLLEAdaptRegisterAllCalled;
73d177a5cSEmil Constantinescu static PetscClassId      TSGLLEADAPT_CLASSID;
83d177a5cSEmil Constantinescu 
93d177a5cSEmil Constantinescu struct _TSGLLEAdaptOps {
103d177a5cSEmil Constantinescu   PetscErrorCode (*choose)(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
113d177a5cSEmil Constantinescu   PetscErrorCode (*destroy)(TSGLLEAdapt);
123d177a5cSEmil Constantinescu   PetscErrorCode (*view)(TSGLLEAdapt, PetscViewer);
13dbbe0bcdSBarry Smith   PetscErrorCode (*setfromoptions)(TSGLLEAdapt, PetscOptionItems *);
143d177a5cSEmil Constantinescu };
153d177a5cSEmil Constantinescu 
163d177a5cSEmil Constantinescu struct _p_TSGLLEAdapt {
173d177a5cSEmil Constantinescu   PETSCHEADER(struct _TSGLLEAdaptOps);
183d177a5cSEmil Constantinescu   void *data;
193d177a5cSEmil Constantinescu };
203d177a5cSEmil Constantinescu 
213d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt);
223d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt);
233d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt);
243d177a5cSEmil Constantinescu 
253d177a5cSEmil Constantinescu /*@C
263d177a5cSEmil Constantinescu    TSGLLEAdaptRegister -  adds a TSGLLEAdapt implementation
273d177a5cSEmil Constantinescu 
283d177a5cSEmil Constantinescu    Not Collective
293d177a5cSEmil Constantinescu 
303d177a5cSEmil Constantinescu    Input Parameters:
313d177a5cSEmil Constantinescu +  name_scheme - name of user-defined adaptivity scheme
323d177a5cSEmil Constantinescu -  routine_create - routine to create method context
333d177a5cSEmil Constantinescu 
343d177a5cSEmil Constantinescu    Notes:
353d177a5cSEmil Constantinescu    TSGLLEAdaptRegister() may be called multiple times to add several user-defined families.
363d177a5cSEmil Constantinescu 
373d177a5cSEmil Constantinescu    Sample usage:
383d177a5cSEmil Constantinescu .vb
393d177a5cSEmil Constantinescu    TSGLLEAdaptRegister("my_scheme",MySchemeCreate);
403d177a5cSEmil Constantinescu .ve
413d177a5cSEmil Constantinescu 
423d177a5cSEmil Constantinescu    Then, your scheme can be chosen with the procedural interface via
433d177a5cSEmil Constantinescu $     TSGLLEAdaptSetType(ts,"my_scheme")
443d177a5cSEmil Constantinescu    or at runtime via the option
453d177a5cSEmil Constantinescu $     -ts_adapt_type my_scheme
463d177a5cSEmil Constantinescu 
473d177a5cSEmil Constantinescu    Level: advanced
483d177a5cSEmil Constantinescu 
49db781477SPatrick Sanan .seealso: `TSGLLEAdaptRegisterAll()`
503d177a5cSEmil Constantinescu @*/
519371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptRegister(const char sname[], PetscErrorCode (*function)(TSGLLEAdapt)) {
523d177a5cSEmil Constantinescu   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptInitializePackage());
549566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSGLLEAdaptList, sname, function));
553d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
563d177a5cSEmil Constantinescu }
573d177a5cSEmil Constantinescu 
583d177a5cSEmil Constantinescu /*@C
593d177a5cSEmil Constantinescu   TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLLEAdapt
603d177a5cSEmil Constantinescu 
613d177a5cSEmil Constantinescu   Not Collective
623d177a5cSEmil Constantinescu 
633d177a5cSEmil Constantinescu   Level: advanced
643d177a5cSEmil Constantinescu 
65db781477SPatrick Sanan .seealso: `TSGLLEAdaptRegisterDestroy()`
663d177a5cSEmil Constantinescu @*/
679371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptRegisterAll(void) {
683d177a5cSEmil Constantinescu   PetscFunctionBegin;
693d177a5cSEmil Constantinescu   if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(0);
703d177a5cSEmil Constantinescu   TSGLLEAdaptRegisterAllCalled = PETSC_TRUE;
719566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_NONE, TSGLLEAdaptCreate_None));
729566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_SIZE, TSGLLEAdaptCreate_Size));
739566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_BOTH, TSGLLEAdaptCreate_Both));
743d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
753d177a5cSEmil Constantinescu }
763d177a5cSEmil Constantinescu 
773d177a5cSEmil Constantinescu /*@C
783d177a5cSEmil Constantinescu   TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is
793d177a5cSEmil Constantinescu   called from PetscFinalize().
803d177a5cSEmil Constantinescu 
813d177a5cSEmil Constantinescu   Level: developer
823d177a5cSEmil Constantinescu 
83db781477SPatrick Sanan .seealso: `PetscFinalize()`
843d177a5cSEmil Constantinescu @*/
859371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptFinalizePackage(void) {
863d177a5cSEmil Constantinescu   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSGLLEAdaptList));
883d177a5cSEmil Constantinescu   TSGLLEAdaptPackageInitialized = PETSC_FALSE;
893d177a5cSEmil Constantinescu   TSGLLEAdaptRegisterAllCalled  = PETSC_FALSE;
903d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
913d177a5cSEmil Constantinescu }
923d177a5cSEmil Constantinescu 
933d177a5cSEmil Constantinescu /*@C
943d177a5cSEmil Constantinescu   TSGLLEAdaptInitializePackage - This function initializes everything in the TSGLLEAdapt package. It is
958a690491SBarry Smith   called from TSInitializePackage().
963d177a5cSEmil Constantinescu 
973d177a5cSEmil Constantinescu   Level: developer
983d177a5cSEmil Constantinescu 
99db781477SPatrick Sanan .seealso: `PetscInitialize()`
1003d177a5cSEmil Constantinescu @*/
1019371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptInitializePackage(void) {
1023d177a5cSEmil Constantinescu   PetscFunctionBegin;
1033d177a5cSEmil Constantinescu   if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(0);
1043d177a5cSEmil Constantinescu   TSGLLEAdaptPackageInitialized = PETSC_TRUE;
1059566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("TSGLLEAdapt", &TSGLLEADAPT_CLASSID));
1069566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegisterAll());
1079566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSGLLEAdaptFinalizePackage));
1083d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1093d177a5cSEmil Constantinescu }
1103d177a5cSEmil Constantinescu 
1119371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt, TSGLLEAdaptType type) {
1125f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TSGLLEAdapt);
1133d177a5cSEmil Constantinescu 
1143d177a5cSEmil Constantinescu   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSGLLEAdaptList, type, &r));
1163c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAdapt type \"%s\" given", type);
117dbbe0bcdSBarry Smith   if (((PetscObject)adapt)->type_name) PetscUseTypeMethod(adapt, destroy);
1189566063dSJacob Faibussowitsch   PetscCall((*r)(adapt));
1199566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type));
1203d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1213d177a5cSEmil Constantinescu }
1223d177a5cSEmil Constantinescu 
1239371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt, const char prefix[]) {
1243d177a5cSEmil Constantinescu   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix));
1263d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1273d177a5cSEmil Constantinescu }
1283d177a5cSEmil Constantinescu 
1299371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt, PetscViewer viewer) {
1303d177a5cSEmil Constantinescu   PetscBool iascii;
1313d177a5cSEmil Constantinescu 
1323d177a5cSEmil Constantinescu   PetscFunctionBegin;
1339566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1343d177a5cSEmil Constantinescu   if (iascii) {
1359566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer));
1363d177a5cSEmil Constantinescu     if (adapt->ops->view) {
1379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
138dbbe0bcdSBarry Smith       PetscUseTypeMethod(adapt, view, viewer);
1399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1403d177a5cSEmil Constantinescu     }
1413d177a5cSEmil Constantinescu   }
1423d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1433d177a5cSEmil Constantinescu }
1443d177a5cSEmil Constantinescu 
1459371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) {
1463d177a5cSEmil Constantinescu   PetscFunctionBegin;
1473d177a5cSEmil Constantinescu   if (!*adapt) PetscFunctionReturn(0);
1483d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(*adapt, TSGLLEADAPT_CLASSID, 1);
1499371c9d4SSatish Balay   if (--((PetscObject)(*adapt))->refct > 0) {
1509371c9d4SSatish Balay     *adapt = NULL;
1519371c9d4SSatish Balay     PetscFunctionReturn(0);
1529371c9d4SSatish Balay   }
153dbbe0bcdSBarry Smith   PetscTryTypeMethod((*adapt), destroy);
1549566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(adapt));
1553d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1563d177a5cSEmil Constantinescu }
1573d177a5cSEmil Constantinescu 
1589371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt adapt, PetscOptionItems *PetscOptionsObject) {
1593d177a5cSEmil Constantinescu   char      type[256] = TSGLLEADAPT_BOTH;
1603d177a5cSEmil Constantinescu   PetscBool flg;
1613d177a5cSEmil Constantinescu 
1623d177a5cSEmil Constantinescu   PetscFunctionBegin;
1633d177a5cSEmil Constantinescu   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this
1643d177a5cSEmil Constantinescu   * function can only be called from inside TSSetFromOptions_GLLE()  */
165d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TSGLLE Adaptivity options");
1669371c9d4SSatish Balay   PetscCall(PetscOptionsFList("-ts_adapt_type", "Algorithm to use for adaptivity", "TSGLLEAdaptSetType", TSGLLEAdaptList, ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type, type, sizeof(type), &flg));
1679566063dSJacob Faibussowitsch   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSGLLEAdaptSetType(adapt, type));
168dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
169d0609cedSBarry Smith   PetscOptionsHeadEnd();
1703d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1713d177a5cSEmil Constantinescu }
1723d177a5cSEmil Constantinescu 
1739371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt adapt, PetscInt n, const PetscInt orders[], const PetscReal errors[], const PetscReal cost[], PetscInt cur, PetscReal h, PetscReal tleft, PetscInt *next_sc, PetscReal *next_h, PetscBool *finish) {
1743d177a5cSEmil Constantinescu   PetscFunctionBegin;
1753d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSGLLEADAPT_CLASSID, 1);
1763d177a5cSEmil Constantinescu   PetscValidIntPointer(orders, 3);
177dadcf809SJacob Faibussowitsch   PetscValidRealPointer(errors, 4);
178dadcf809SJacob Faibussowitsch   PetscValidRealPointer(cost, 5);
1793d177a5cSEmil Constantinescu   PetscValidIntPointer(next_sc, 9);
180dadcf809SJacob Faibussowitsch   PetscValidRealPointer(next_h, 10);
181064a246eSJacob Faibussowitsch   PetscValidBoolPointer(finish, 11);
182dbbe0bcdSBarry Smith   PetscUseTypeMethod(adapt, choose, n, orders, errors, cost, cur, h, tleft, next_sc, next_h, finish);
1833d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1843d177a5cSEmil Constantinescu }
1853d177a5cSEmil Constantinescu 
1869371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm, TSGLLEAdapt *inadapt) {
1873d177a5cSEmil Constantinescu   TSGLLEAdapt adapt;
1883d177a5cSEmil Constantinescu 
1893d177a5cSEmil Constantinescu   PetscFunctionBegin;
1903d177a5cSEmil Constantinescu   *inadapt = NULL;
1919566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(adapt, TSGLLEADAPT_CLASSID, "TSGLLEAdapt", "General Linear adaptivity", "TS", comm, TSGLLEAdaptDestroy, TSGLLEAdaptView));
1923d177a5cSEmil Constantinescu   *inadapt = adapt;
1933d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1943d177a5cSEmil Constantinescu }
1953d177a5cSEmil Constantinescu 
1963d177a5cSEmil Constantinescu /*
1970e3d61c9SBarry Smith    Implementations
1983d177a5cSEmil Constantinescu */
1993d177a5cSEmil Constantinescu 
2009371c9d4SSatish Balay static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) {
2013d177a5cSEmil Constantinescu   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(adapt->data));
2033d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2043d177a5cSEmil Constantinescu }
2053d177a5cSEmil Constantinescu 
2063d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */
2073d177a5cSEmil Constantinescu typedef struct {
2083d177a5cSEmil Constantinescu   PetscInt  scheme;
2093d177a5cSEmil Constantinescu   PetscReal h;
2103d177a5cSEmil Constantinescu } TSGLLEAdapt_None;
2113d177a5cSEmil Constantinescu 
2129371c9d4SSatish Balay static PetscErrorCode TSGLLEAdaptChoose_None(TSGLLEAdapt adapt, PetscInt n, const PetscInt orders[], const PetscReal errors[], const PetscReal cost[], PetscInt cur, PetscReal h, PetscReal tleft, PetscInt *next_sc, PetscReal *next_h, PetscBool *finish) {
2133d177a5cSEmil Constantinescu   PetscFunctionBegin;
2143d177a5cSEmil Constantinescu   *next_sc = cur;
2153d177a5cSEmil Constantinescu   *next_h  = h;
2163d177a5cSEmil Constantinescu   if (*next_h > tleft) {
2173d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
2183d177a5cSEmil Constantinescu     *next_h = tleft;
2193d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
2203d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2213d177a5cSEmil Constantinescu }
2223d177a5cSEmil Constantinescu 
2239371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) {
2243d177a5cSEmil Constantinescu   TSGLLEAdapt_None *a;
2253d177a5cSEmil Constantinescu 
2263d177a5cSEmil Constantinescu   PetscFunctionBegin;
227*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
2283d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
2293d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_None;
2303d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
2313d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2323d177a5cSEmil Constantinescu }
2333d177a5cSEmil Constantinescu 
2343d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */
2353d177a5cSEmil Constantinescu typedef struct {
2363d177a5cSEmil Constantinescu   PetscReal desired_h;
2373d177a5cSEmil Constantinescu } TSGLLEAdapt_Size;
2383d177a5cSEmil Constantinescu 
2399371c9d4SSatish Balay static PetscErrorCode TSGLLEAdaptChoose_Size(TSGLLEAdapt adapt, PetscInt n, const PetscInt orders[], const PetscReal errors[], const PetscReal cost[], PetscInt cur, PetscReal h, PetscReal tleft, PetscInt *next_sc, PetscReal *next_h, PetscBool *finish) {
2403d177a5cSEmil Constantinescu   TSGLLEAdapt_Size *sz  = (TSGLLEAdapt_Size *)adapt->data;
2413d177a5cSEmil Constantinescu   PetscReal         dec = 0.2, inc = 5.0, safe = 0.9, optimal, last_desired_h;
2423d177a5cSEmil Constantinescu 
2433d177a5cSEmil Constantinescu   PetscFunctionBegin;
2443d177a5cSEmil Constantinescu   *next_sc       = cur;
2453d177a5cSEmil Constantinescu   optimal        = PetscPowReal((PetscReal)errors[cur], (PetscReal)-1. / (safe * orders[cur]));
2463d177a5cSEmil Constantinescu   /* Step sizes oscillate when there is no smoothing.  Here we use a geometric mean of the current step size and the
2473d177a5cSEmil Constantinescu   * one that would have been taken (without smoothing) on the last step. */
2483d177a5cSEmil Constantinescu   last_desired_h = sz->desired_h;
2493d177a5cSEmil Constantinescu   sz->desired_h  = h * PetscMax(dec, PetscMin(inc, optimal)); /* Trim to [dec,inc] */
2503d177a5cSEmil Constantinescu 
2513d177a5cSEmil Constantinescu   /* Normally only happens on the first step */
2523d177a5cSEmil Constantinescu   if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h);
2533d177a5cSEmil Constantinescu   else *next_h = sz->desired_h;
2543d177a5cSEmil Constantinescu 
2553d177a5cSEmil Constantinescu   if (*next_h > tleft) {
2563d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
2573d177a5cSEmil Constantinescu     *next_h = tleft;
2583d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
2593d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2603d177a5cSEmil Constantinescu }
2613d177a5cSEmil Constantinescu 
2629371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) {
2633d177a5cSEmil Constantinescu   TSGLLEAdapt_Size *a;
2643d177a5cSEmil Constantinescu 
2653d177a5cSEmil Constantinescu   PetscFunctionBegin;
266*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
2673d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
2683d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_Size;
2693d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
2703d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2713d177a5cSEmil Constantinescu }
2723d177a5cSEmil Constantinescu 
2733d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */
2743d177a5cSEmil Constantinescu typedef struct {
2753d177a5cSEmil Constantinescu   PetscInt  count_at_order;
2763d177a5cSEmil Constantinescu   PetscReal desired_h;
2773d177a5cSEmil Constantinescu } TSGLLEAdapt_Both;
2783d177a5cSEmil Constantinescu 
2799371c9d4SSatish Balay static PetscErrorCode TSGLLEAdaptChoose_Both(TSGLLEAdapt adapt, PetscInt n, const PetscInt orders[], const PetscReal errors[], const PetscReal cost[], PetscInt cur, PetscReal h, PetscReal tleft, PetscInt *next_sc, PetscReal *next_h, PetscBool *finish) {
2803d177a5cSEmil Constantinescu   TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both *)adapt->data;
2813d177a5cSEmil Constantinescu   PetscReal         dec = 0.2, inc = 5.0, safe = 0.9;
2829371c9d4SSatish Balay   struct {
2839371c9d4SSatish Balay     PetscInt  id;
2849371c9d4SSatish Balay     PetscReal h, eff;
2859371c9d4SSatish Balay   } best = {-1, 0, 0}, trial = {-1, 0, 0}, current = {-1, 0, 0};
2863d177a5cSEmil Constantinescu   PetscInt i;
2873d177a5cSEmil Constantinescu 
2883d177a5cSEmil Constantinescu   PetscFunctionBegin;
2893d177a5cSEmil Constantinescu   for (i = 0; i < n; i++) {
2903d177a5cSEmil Constantinescu     PetscReal optimal;
2913d177a5cSEmil Constantinescu     trial.id  = i;
2923d177a5cSEmil Constantinescu     optimal   = PetscPowReal((PetscReal)errors[i], (PetscReal)-1. / (safe * orders[i]));
2933d177a5cSEmil Constantinescu     trial.h   = h * optimal;
2943d177a5cSEmil Constantinescu     trial.eff = trial.h / cost[i];
2959566063dSJacob Faibussowitsch     if (trial.eff > best.eff) PetscCall(PetscArraycpy(&best, &trial, 1));
2969566063dSJacob Faibussowitsch     if (i == cur) PetscCall(PetscArraycpy(&current, &trial, 1));
2973d177a5cSEmil Constantinescu   }
2983d177a5cSEmil Constantinescu   /* Only switch orders if the scheme offers significant benefits over the current one.
2993d177a5cSEmil Constantinescu   When the scheme is not changing, only change step size if it offers significant benefits. */
3003d177a5cSEmil Constantinescu   if (best.eff < 1.2 * current.eff || both->count_at_order < orders[cur] + 2) {
3013d177a5cSEmil Constantinescu     PetscReal last_desired_h;
3023d177a5cSEmil Constantinescu     *next_sc        = current.id;
3033d177a5cSEmil Constantinescu     last_desired_h  = both->desired_h;
3043d177a5cSEmil Constantinescu     both->desired_h = PetscMax(h * dec, PetscMin(h * inc, current.h));
3059371c9d4SSatish Balay     *next_h         = (both->count_at_order > 0) ? PetscSqrtReal(last_desired_h * both->desired_h) : both->desired_h;
3063d177a5cSEmil Constantinescu     both->count_at_order++;
3073d177a5cSEmil Constantinescu   } else {
3083d177a5cSEmil Constantinescu     PetscReal rat        = cost[best.id] / cost[cur];
3093d177a5cSEmil Constantinescu     *next_sc             = best.id;
3103d177a5cSEmil Constantinescu     *next_h              = PetscMax(h * rat * dec, PetscMin(h * rat * inc, best.h));
3113d177a5cSEmil Constantinescu     both->count_at_order = 0;
3123d177a5cSEmil Constantinescu     both->desired_h      = best.h;
3133d177a5cSEmil Constantinescu   }
3143d177a5cSEmil Constantinescu 
3153d177a5cSEmil Constantinescu   if (*next_h > tleft) {
3163d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
3173d177a5cSEmil Constantinescu     *next_h = tleft;
3183d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
3193d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3203d177a5cSEmil Constantinescu }
3213d177a5cSEmil Constantinescu 
3229371c9d4SSatish Balay PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) {
3233d177a5cSEmil Constantinescu   TSGLLEAdapt_Both *a;
3243d177a5cSEmil Constantinescu 
3253d177a5cSEmil Constantinescu   PetscFunctionBegin;
326*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
3273d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
3283d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_Both;
3293d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
3303d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3313d177a5cSEmil Constantinescu }
332