xref: /petsc/src/ts/impls/implicit/glle/glleadapt.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
13d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I  "petscts.h" I*/
23d177a5cSEmil Constantinescu 
33d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAdaptList;
43d177a5cSEmil Constantinescu static PetscBool         TSGLLEAdaptPackageInitialized;
53d177a5cSEmil Constantinescu static PetscBool         TSGLLEAdaptRegisterAllCalled;
63d177a5cSEmil Constantinescu static PetscClassId      TSGLLEADAPT_CLASSID;
73d177a5cSEmil Constantinescu 
83d177a5cSEmil Constantinescu struct _TSGLLEAdaptOps {
93d177a5cSEmil Constantinescu   PetscErrorCode (*choose)(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *);
103d177a5cSEmil Constantinescu   PetscErrorCode (*destroy)(TSGLLEAdapt);
113d177a5cSEmil Constantinescu   PetscErrorCode (*view)(TSGLLEAdapt, PetscViewer);
12dbbe0bcdSBarry Smith   PetscErrorCode (*setfromoptions)(TSGLLEAdapt, PetscOptionItems *);
133d177a5cSEmil Constantinescu };
143d177a5cSEmil Constantinescu 
153d177a5cSEmil Constantinescu struct _p_TSGLLEAdapt {
163d177a5cSEmil Constantinescu   PETSCHEADER(struct _TSGLLEAdaptOps);
173d177a5cSEmil Constantinescu   void *data;
183d177a5cSEmil Constantinescu };
193d177a5cSEmil Constantinescu 
203d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt);
213d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt);
223d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt);
233d177a5cSEmil Constantinescu 
243d177a5cSEmil Constantinescu /*@C
25bcf0153eSBarry Smith   TSGLLEAdaptRegister -  adds a `TSGLLEAdapt` implementation
263d177a5cSEmil Constantinescu 
273d177a5cSEmil Constantinescu   Not Collective
283d177a5cSEmil Constantinescu 
293d177a5cSEmil Constantinescu   Input Parameters:
3020f4b53cSBarry Smith + sname    - name of user-defined adaptivity scheme
3120f4b53cSBarry Smith - function - routine to create method context
3220f4b53cSBarry Smith 
3320f4b53cSBarry Smith   Level: advanced
343d177a5cSEmil Constantinescu 
35bcf0153eSBarry Smith   Note:
36bcf0153eSBarry Smith   `TSGLLEAdaptRegister()` may be called multiple times to add several user-defined families.
373d177a5cSEmil Constantinescu 
38b43aa488SJacob Faibussowitsch   Example Usage:
393d177a5cSEmil Constantinescu .vb
403d177a5cSEmil Constantinescu    TSGLLEAdaptRegister("my_scheme", MySchemeCreate);
413d177a5cSEmil Constantinescu .ve
423d177a5cSEmil Constantinescu 
433d177a5cSEmil Constantinescu   Then, your scheme can be chosen with the procedural interface via
443d177a5cSEmil Constantinescu $     TSGLLEAdaptSetType(ts, "my_scheme")
453d177a5cSEmil Constantinescu   or at runtime via the option
463d177a5cSEmil Constantinescu $     -ts_adapt_type my_scheme
473d177a5cSEmil Constantinescu 
481cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLEAdapt`, `TSGLLEAdaptRegisterAll()`
493d177a5cSEmil Constantinescu @*/
50d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegister(const char sname[], PetscErrorCode (*function)(TSGLLEAdapt))
51d71ae5a4SJacob Faibussowitsch {
523d177a5cSEmil Constantinescu   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptInitializePackage());
549566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSGLLEAdaptList, sname, function));
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563d177a5cSEmil Constantinescu }
573d177a5cSEmil Constantinescu 
583d177a5cSEmil Constantinescu /*@C
59bcf0153eSBarry Smith   TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in `TSGLLEAdapt`
603d177a5cSEmil Constantinescu 
613d177a5cSEmil Constantinescu   Not Collective
623d177a5cSEmil Constantinescu 
633d177a5cSEmil Constantinescu   Level: advanced
643d177a5cSEmil Constantinescu 
651cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdapt`, `TSGLLE`, `TSGLLEAdaptRegisterDestroy()`
663d177a5cSEmil Constantinescu @*/
67d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegisterAll(void)
68d71ae5a4SJacob Faibussowitsch {
693d177a5cSEmil Constantinescu   PetscFunctionBegin;
703ba16761SJacob Faibussowitsch   if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
713d177a5cSEmil Constantinescu   TSGLLEAdaptRegisterAllCalled = PETSC_TRUE;
729566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_NONE, TSGLLEAdaptCreate_None));
739566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_SIZE, TSGLLEAdaptCreate_Size));
749566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_BOTH, TSGLLEAdaptCreate_Both));
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
763d177a5cSEmil Constantinescu }
773d177a5cSEmil Constantinescu 
783d177a5cSEmil Constantinescu /*@C
79b43aa488SJacob Faibussowitsch   TSGLLEAdaptFinalizePackage - This function destroys everything in the `TSGLLE` package. It is
80bcf0153eSBarry Smith   called from `PetscFinalize()`.
813d177a5cSEmil Constantinescu 
823d177a5cSEmil Constantinescu   Level: developer
833d177a5cSEmil Constantinescu 
841cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSGLLEAdapt`, `TSGLLEAdaptInitializePackage()`
853d177a5cSEmil Constantinescu @*/
86d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptFinalizePackage(void)
87d71ae5a4SJacob Faibussowitsch {
883d177a5cSEmil Constantinescu   PetscFunctionBegin;
899566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSGLLEAdaptList));
903d177a5cSEmil Constantinescu   TSGLLEAdaptPackageInitialized = PETSC_FALSE;
913d177a5cSEmil Constantinescu   TSGLLEAdaptRegisterAllCalled  = PETSC_FALSE;
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
933d177a5cSEmil Constantinescu }
943d177a5cSEmil Constantinescu 
953d177a5cSEmil Constantinescu /*@C
96bcf0153eSBarry Smith   TSGLLEAdaptInitializePackage - This function initializes everything in the `TSGLLEAdapt` package. It is
97bcf0153eSBarry Smith   called from `TSInitializePackage()`.
983d177a5cSEmil Constantinescu 
993d177a5cSEmil Constantinescu   Level: developer
1003d177a5cSEmil Constantinescu 
1011cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSGLLEAdapt`, `TSGLLEAdaptFinalizePackage()`
1023d177a5cSEmil Constantinescu @*/
103d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptInitializePackage(void)
104d71ae5a4SJacob Faibussowitsch {
1053d177a5cSEmil Constantinescu   PetscFunctionBegin;
1063ba16761SJacob Faibussowitsch   if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1073d177a5cSEmil Constantinescu   TSGLLEAdaptPackageInitialized = PETSC_TRUE;
1089566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("TSGLLEAdapt", &TSGLLEADAPT_CLASSID));
1099566063dSJacob Faibussowitsch   PetscCall(TSGLLEAdaptRegisterAll());
1109566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSGLLEAdaptFinalizePackage));
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1123d177a5cSEmil Constantinescu }
1133d177a5cSEmil Constantinescu 
114d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt, TSGLLEAdaptType type)
115d71ae5a4SJacob Faibussowitsch {
1165f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TSGLLEAdapt);
1173d177a5cSEmil Constantinescu 
1183d177a5cSEmil Constantinescu   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSGLLEAdaptList, type, &r));
1203c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAdapt type \"%s\" given", type);
121dbbe0bcdSBarry Smith   if (((PetscObject)adapt)->type_name) PetscUseTypeMethod(adapt, destroy);
1229566063dSJacob Faibussowitsch   PetscCall((*r)(adapt));
1239566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253d177a5cSEmil Constantinescu }
1263d177a5cSEmil Constantinescu 
127d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt, const char prefix[])
128d71ae5a4SJacob Faibussowitsch {
1293d177a5cSEmil Constantinescu   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix));
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1323d177a5cSEmil Constantinescu }
1333d177a5cSEmil Constantinescu 
134d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt, PetscViewer viewer)
135d71ae5a4SJacob Faibussowitsch {
1363d177a5cSEmil Constantinescu   PetscBool iascii;
1373d177a5cSEmil Constantinescu 
1383d177a5cSEmil Constantinescu   PetscFunctionBegin;
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1403d177a5cSEmil Constantinescu   if (iascii) {
1419566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer));
1423d177a5cSEmil Constantinescu     if (adapt->ops->view) {
1439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
144dbbe0bcdSBarry Smith       PetscUseTypeMethod(adapt, view, viewer);
1459566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1463d177a5cSEmil Constantinescu     }
1473d177a5cSEmil Constantinescu   }
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1493d177a5cSEmil Constantinescu }
1503d177a5cSEmil Constantinescu 
151d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt)
152d71ae5a4SJacob Faibussowitsch {
1533d177a5cSEmil Constantinescu   PetscFunctionBegin;
1543ba16761SJacob Faibussowitsch   if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS);
1553d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(*adapt, TSGLLEADAPT_CLASSID, 1);
156*f4f49eeaSPierre Jolivet   if (--((PetscObject)*adapt)->refct > 0) {
1579371c9d4SSatish Balay     *adapt = NULL;
1583ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1599371c9d4SSatish Balay   }
160*f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*adapt, destroy);
1619566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(adapt));
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1633d177a5cSEmil Constantinescu }
1643d177a5cSEmil Constantinescu 
165d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt adapt, PetscOptionItems *PetscOptionsObject)
166d71ae5a4SJacob Faibussowitsch {
1673d177a5cSEmil Constantinescu   char      type[256] = TSGLLEADAPT_BOTH;
1683d177a5cSEmil Constantinescu   PetscBool flg;
1693d177a5cSEmil Constantinescu 
1703d177a5cSEmil Constantinescu   PetscFunctionBegin;
1713d177a5cSEmil Constantinescu   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this
1723d177a5cSEmil Constantinescu   * function can only be called from inside TSSetFromOptions_GLLE()  */
173d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TSGLLE Adaptivity options");
1749371c9d4SSatish 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));
1759566063dSJacob Faibussowitsch   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSGLLEAdaptSetType(adapt, type));
176dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
177d0609cedSBarry Smith   PetscOptionsHeadEnd();
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1793d177a5cSEmil Constantinescu }
1803d177a5cSEmil Constantinescu 
181d71ae5a4SJacob Faibussowitsch 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)
182d71ae5a4SJacob Faibussowitsch {
1833d177a5cSEmil Constantinescu   PetscFunctionBegin;
1843d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSGLLEADAPT_CLASSID, 1);
1854f572ea9SToby Isaac   PetscAssertPointer(orders, 3);
1864f572ea9SToby Isaac   PetscAssertPointer(errors, 4);
1874f572ea9SToby Isaac   PetscAssertPointer(cost, 5);
1884f572ea9SToby Isaac   PetscAssertPointer(next_sc, 9);
1894f572ea9SToby Isaac   PetscAssertPointer(next_h, 10);
1904f572ea9SToby Isaac   PetscAssertPointer(finish, 11);
191dbbe0bcdSBarry Smith   PetscUseTypeMethod(adapt, choose, n, orders, errors, cost, cur, h, tleft, next_sc, next_h, finish);
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1933d177a5cSEmil Constantinescu }
1943d177a5cSEmil Constantinescu 
195d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm, TSGLLEAdapt *inadapt)
196d71ae5a4SJacob Faibussowitsch {
1973d177a5cSEmil Constantinescu   TSGLLEAdapt adapt;
1983d177a5cSEmil Constantinescu 
1993d177a5cSEmil Constantinescu   PetscFunctionBegin;
2003d177a5cSEmil Constantinescu   *inadapt = NULL;
2019566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(adapt, TSGLLEADAPT_CLASSID, "TSGLLEAdapt", "General Linear adaptivity", "TS", comm, TSGLLEAdaptDestroy, TSGLLEAdaptView));
2023d177a5cSEmil Constantinescu   *inadapt = adapt;
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2043d177a5cSEmil Constantinescu }
2053d177a5cSEmil Constantinescu 
2063d177a5cSEmil Constantinescu /*
2070e3d61c9SBarry Smith    Implementations
2083d177a5cSEmil Constantinescu */
2093d177a5cSEmil Constantinescu 
210d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt)
211d71ae5a4SJacob Faibussowitsch {
2123d177a5cSEmil Constantinescu   PetscFunctionBegin;
2139566063dSJacob Faibussowitsch   PetscCall(PetscFree(adapt->data));
2143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2153d177a5cSEmil Constantinescu }
2163d177a5cSEmil Constantinescu 
2173d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */
2183d177a5cSEmil Constantinescu typedef struct {
2193d177a5cSEmil Constantinescu   PetscInt  scheme;
2203d177a5cSEmil Constantinescu   PetscReal h;
2213d177a5cSEmil Constantinescu } TSGLLEAdapt_None;
2223d177a5cSEmil Constantinescu 
223d71ae5a4SJacob Faibussowitsch 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)
224d71ae5a4SJacob Faibussowitsch {
2253d177a5cSEmil Constantinescu   PetscFunctionBegin;
2263d177a5cSEmil Constantinescu   *next_sc = cur;
2273d177a5cSEmil Constantinescu   *next_h  = h;
2283d177a5cSEmil Constantinescu   if (*next_h > tleft) {
2293d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
2303d177a5cSEmil Constantinescu     *next_h = tleft;
2313d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2333d177a5cSEmil Constantinescu }
2343d177a5cSEmil Constantinescu 
235d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt)
236d71ae5a4SJacob Faibussowitsch {
2373d177a5cSEmil Constantinescu   TSGLLEAdapt_None *a;
2383d177a5cSEmil Constantinescu 
2393d177a5cSEmil Constantinescu   PetscFunctionBegin;
2404dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
2413d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
2423d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_None;
2433d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2453d177a5cSEmil Constantinescu }
2463d177a5cSEmil Constantinescu 
2473d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */
2483d177a5cSEmil Constantinescu typedef struct {
2493d177a5cSEmil Constantinescu   PetscReal desired_h;
2503d177a5cSEmil Constantinescu } TSGLLEAdapt_Size;
2513d177a5cSEmil Constantinescu 
252d71ae5a4SJacob Faibussowitsch 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)
253d71ae5a4SJacob Faibussowitsch {
2543d177a5cSEmil Constantinescu   TSGLLEAdapt_Size *sz  = (TSGLLEAdapt_Size *)adapt->data;
2553d177a5cSEmil Constantinescu   PetscReal         dec = 0.2, inc = 5.0, safe = 0.9, optimal, last_desired_h;
2563d177a5cSEmil Constantinescu 
2573d177a5cSEmil Constantinescu   PetscFunctionBegin;
2583d177a5cSEmil Constantinescu   *next_sc = cur;
2593d177a5cSEmil Constantinescu   optimal  = PetscPowReal((PetscReal)errors[cur], (PetscReal)-1. / (safe * orders[cur]));
2603d177a5cSEmil Constantinescu   /* Step sizes oscillate when there is no smoothing.  Here we use a geometric mean of the current step size and the
2613d177a5cSEmil Constantinescu   * one that would have been taken (without smoothing) on the last step. */
2623d177a5cSEmil Constantinescu   last_desired_h = sz->desired_h;
2633d177a5cSEmil Constantinescu   sz->desired_h  = h * PetscMax(dec, PetscMin(inc, optimal)); /* Trim to [dec,inc] */
2643d177a5cSEmil Constantinescu 
2653d177a5cSEmil Constantinescu   /* Normally only happens on the first step */
2663d177a5cSEmil Constantinescu   if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h);
2673d177a5cSEmil Constantinescu   else *next_h = sz->desired_h;
2683d177a5cSEmil Constantinescu 
2693d177a5cSEmil Constantinescu   if (*next_h > tleft) {
2703d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
2713d177a5cSEmil Constantinescu     *next_h = tleft;
2723d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2743d177a5cSEmil Constantinescu }
2753d177a5cSEmil Constantinescu 
276d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt)
277d71ae5a4SJacob Faibussowitsch {
2783d177a5cSEmil Constantinescu   TSGLLEAdapt_Size *a;
2793d177a5cSEmil Constantinescu 
2803d177a5cSEmil Constantinescu   PetscFunctionBegin;
2814dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
2823d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
2833d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_Size;
2843d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2863d177a5cSEmil Constantinescu }
2873d177a5cSEmil Constantinescu 
2883d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */
2893d177a5cSEmil Constantinescu typedef struct {
2903d177a5cSEmil Constantinescu   PetscInt  count_at_order;
2913d177a5cSEmil Constantinescu   PetscReal desired_h;
2923d177a5cSEmil Constantinescu } TSGLLEAdapt_Both;
2933d177a5cSEmil Constantinescu 
294d71ae5a4SJacob Faibussowitsch 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)
295d71ae5a4SJacob Faibussowitsch {
2963d177a5cSEmil Constantinescu   TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both *)adapt->data;
2973d177a5cSEmil Constantinescu   PetscReal         dec = 0.2, inc = 5.0, safe = 0.9;
2989371c9d4SSatish Balay   struct {
2999371c9d4SSatish Balay     PetscInt  id;
3009371c9d4SSatish Balay     PetscReal h, eff;
3019371c9d4SSatish Balay   } best = {-1, 0, 0}, trial = {-1, 0, 0}, current = {-1, 0, 0};
3023d177a5cSEmil Constantinescu   PetscInt i;
3033d177a5cSEmil Constantinescu 
3043d177a5cSEmil Constantinescu   PetscFunctionBegin;
3053d177a5cSEmil Constantinescu   for (i = 0; i < n; i++) {
3063d177a5cSEmil Constantinescu     PetscReal optimal;
3073d177a5cSEmil Constantinescu     trial.id  = i;
3083d177a5cSEmil Constantinescu     optimal   = PetscPowReal((PetscReal)errors[i], (PetscReal)-1. / (safe * orders[i]));
3093d177a5cSEmil Constantinescu     trial.h   = h * optimal;
3103d177a5cSEmil Constantinescu     trial.eff = trial.h / cost[i];
3119566063dSJacob Faibussowitsch     if (trial.eff > best.eff) PetscCall(PetscArraycpy(&best, &trial, 1));
3129566063dSJacob Faibussowitsch     if (i == cur) PetscCall(PetscArraycpy(&current, &trial, 1));
3133d177a5cSEmil Constantinescu   }
3143d177a5cSEmil Constantinescu   /* Only switch orders if the scheme offers significant benefits over the current one.
3153d177a5cSEmil Constantinescu   When the scheme is not changing, only change step size if it offers significant benefits. */
3163d177a5cSEmil Constantinescu   if (best.eff < 1.2 * current.eff || both->count_at_order < orders[cur] + 2) {
3173d177a5cSEmil Constantinescu     PetscReal last_desired_h;
3183d177a5cSEmil Constantinescu     *next_sc        = current.id;
3193d177a5cSEmil Constantinescu     last_desired_h  = both->desired_h;
3203d177a5cSEmil Constantinescu     both->desired_h = PetscMax(h * dec, PetscMin(h * inc, current.h));
3219371c9d4SSatish Balay     *next_h         = (both->count_at_order > 0) ? PetscSqrtReal(last_desired_h * both->desired_h) : both->desired_h;
3223d177a5cSEmil Constantinescu     both->count_at_order++;
3233d177a5cSEmil Constantinescu   } else {
3243d177a5cSEmil Constantinescu     PetscReal rat        = cost[best.id] / cost[cur];
3253d177a5cSEmil Constantinescu     *next_sc             = best.id;
3263d177a5cSEmil Constantinescu     *next_h              = PetscMax(h * rat * dec, PetscMin(h * rat * inc, best.h));
3273d177a5cSEmil Constantinescu     both->count_at_order = 0;
3283d177a5cSEmil Constantinescu     both->desired_h      = best.h;
3293d177a5cSEmil Constantinescu   }
3303d177a5cSEmil Constantinescu 
3313d177a5cSEmil Constantinescu   if (*next_h > tleft) {
3323d177a5cSEmil Constantinescu     *finish = PETSC_TRUE;
3333d177a5cSEmil Constantinescu     *next_h = tleft;
3343d177a5cSEmil Constantinescu   } else *finish = PETSC_FALSE;
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3363d177a5cSEmil Constantinescu }
3373d177a5cSEmil Constantinescu 
338d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt)
339d71ae5a4SJacob Faibussowitsch {
3403d177a5cSEmil Constantinescu   TSGLLEAdapt_Both *a;
3413d177a5cSEmil Constantinescu 
3423d177a5cSEmil Constantinescu   PetscFunctionBegin;
3434dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
3443d177a5cSEmil Constantinescu   adapt->data         = (void *)a;
3453d177a5cSEmil Constantinescu   adapt->ops->choose  = TSGLLEAdaptChoose_Both;
3463d177a5cSEmil Constantinescu   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3483d177a5cSEmil Constantinescu }
349