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); 12ce78bad3SBarry 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 27cc4c1da9SBarry Smith Not Collective, No Fortran Support 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 44*b44f4de4SBarry Smith .vb 45*b44f4de4SBarry Smith TSGLLEAdaptSetType(ts, "my_scheme") 46*b44f4de4SBarry Smith .ve 473d177a5cSEmil Constantinescu or at runtime via the option 48*b44f4de4SBarry Smith .vb 49*b44f4de4SBarry Smith -ts_adapt_type my_scheme 50*b44f4de4SBarry Smith .ve 513d177a5cSEmil Constantinescu 521cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLEAdapt`, `TSGLLEAdaptRegisterAll()` 533d177a5cSEmil Constantinescu @*/ 54d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegister(const char sname[], PetscErrorCode (*function)(TSGLLEAdapt)) 55d71ae5a4SJacob Faibussowitsch { 563d177a5cSEmil Constantinescu PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptInitializePackage()); 589566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEAdaptList, sname, function)); 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 603d177a5cSEmil Constantinescu } 613d177a5cSEmil Constantinescu 623d177a5cSEmil Constantinescu /*@C 63bcf0153eSBarry Smith TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in `TSGLLEAdapt` 643d177a5cSEmil Constantinescu 653d177a5cSEmil Constantinescu Not Collective 663d177a5cSEmil Constantinescu 673d177a5cSEmil Constantinescu Level: advanced 683d177a5cSEmil Constantinescu 691cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdapt`, `TSGLLE`, `TSGLLEAdaptRegisterDestroy()` 703d177a5cSEmil Constantinescu @*/ 71d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegisterAll(void) 72d71ae5a4SJacob Faibussowitsch { 733d177a5cSEmil Constantinescu PetscFunctionBegin; 743ba16761SJacob Faibussowitsch if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 753d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_TRUE; 769566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_NONE, TSGLLEAdaptCreate_None)); 779566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_SIZE, TSGLLEAdaptCreate_Size)); 789566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_BOTH, TSGLLEAdaptCreate_Both)); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 803d177a5cSEmil Constantinescu } 813d177a5cSEmil Constantinescu 823d177a5cSEmil Constantinescu /*@C 83b43aa488SJacob Faibussowitsch TSGLLEAdaptFinalizePackage - This function destroys everything in the `TSGLLE` package. It is 84bcf0153eSBarry Smith called from `PetscFinalize()`. 853d177a5cSEmil Constantinescu 863d177a5cSEmil Constantinescu Level: developer 873d177a5cSEmil Constantinescu 881cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSGLLEAdapt`, `TSGLLEAdaptInitializePackage()` 893d177a5cSEmil Constantinescu @*/ 90d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptFinalizePackage(void) 91d71ae5a4SJacob Faibussowitsch { 923d177a5cSEmil Constantinescu PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEAdaptList)); 943d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_FALSE; 953d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_FALSE; 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 973d177a5cSEmil Constantinescu } 983d177a5cSEmil Constantinescu 993d177a5cSEmil Constantinescu /*@C 100bcf0153eSBarry Smith TSGLLEAdaptInitializePackage - This function initializes everything in the `TSGLLEAdapt` package. It is 101bcf0153eSBarry Smith called from `TSInitializePackage()`. 1023d177a5cSEmil Constantinescu 1033d177a5cSEmil Constantinescu Level: developer 1043d177a5cSEmil Constantinescu 1051cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSGLLEAdapt`, `TSGLLEAdaptFinalizePackage()` 1063d177a5cSEmil Constantinescu @*/ 107d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptInitializePackage(void) 108d71ae5a4SJacob Faibussowitsch { 1093d177a5cSEmil Constantinescu PetscFunctionBegin; 1103ba16761SJacob Faibussowitsch if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 1113d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_TRUE; 1129566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSGLLEAdapt", &TSGLLEADAPT_CLASSID)); 1139566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegisterAll()); 1149566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSGLLEAdaptFinalizePackage)); 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1163d177a5cSEmil Constantinescu } 1173d177a5cSEmil Constantinescu 118d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt, TSGLLEAdaptType type) 119d71ae5a4SJacob Faibussowitsch { 1205f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSGLLEAdapt); 1213d177a5cSEmil Constantinescu 1223d177a5cSEmil Constantinescu PetscFunctionBegin; 1239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEAdaptList, type, &r)); 1243c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAdapt type \"%s\" given", type); 125dbbe0bcdSBarry Smith if (((PetscObject)adapt)->type_name) PetscUseTypeMethod(adapt, destroy); 1269566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 1279566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 1283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1293d177a5cSEmil Constantinescu } 1303d177a5cSEmil Constantinescu 131d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt, const char prefix[]) 132d71ae5a4SJacob Faibussowitsch { 1333d177a5cSEmil Constantinescu PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1363d177a5cSEmil Constantinescu } 1373d177a5cSEmil Constantinescu 138d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt, PetscViewer viewer) 139d71ae5a4SJacob Faibussowitsch { 1403d177a5cSEmil Constantinescu PetscBool iascii; 1413d177a5cSEmil Constantinescu 1423d177a5cSEmil Constantinescu PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1443d177a5cSEmil Constantinescu if (iascii) { 1459566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 1463d177a5cSEmil Constantinescu if (adapt->ops->view) { 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 148dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, view, viewer); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1503d177a5cSEmil Constantinescu } 1513d177a5cSEmil Constantinescu } 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1533d177a5cSEmil Constantinescu } 1543d177a5cSEmil Constantinescu 155d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) 156d71ae5a4SJacob Faibussowitsch { 1573d177a5cSEmil Constantinescu PetscFunctionBegin; 1583ba16761SJacob Faibussowitsch if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS); 1593d177a5cSEmil Constantinescu PetscValidHeaderSpecific(*adapt, TSGLLEADAPT_CLASSID, 1); 160f4f49eeaSPierre Jolivet if (--((PetscObject)*adapt)->refct > 0) { 1619371c9d4SSatish Balay *adapt = NULL; 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1639371c9d4SSatish Balay } 164f4f49eeaSPierre Jolivet PetscTryTypeMethod(*adapt, destroy); 1659566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1673d177a5cSEmil Constantinescu } 1683d177a5cSEmil Constantinescu 169ce78bad3SBarry Smith PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt adapt, PetscOptionItems PetscOptionsObject) 170d71ae5a4SJacob Faibussowitsch { 1713d177a5cSEmil Constantinescu char type[256] = TSGLLEADAPT_BOTH; 1723d177a5cSEmil Constantinescu PetscBool flg; 1733d177a5cSEmil Constantinescu 1743d177a5cSEmil Constantinescu PetscFunctionBegin; 1753d177a5cSEmil Constantinescu /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this 1763d177a5cSEmil Constantinescu * function can only be called from inside TSSetFromOptions_GLLE() */ 177d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TSGLLE Adaptivity options"); 1789371c9d4SSatish 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)); 1799566063dSJacob Faibussowitsch if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSGLLEAdaptSetType(adapt, type)); 180dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 181d0609cedSBarry Smith PetscOptionsHeadEnd(); 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1833d177a5cSEmil Constantinescu } 1843d177a5cSEmil Constantinescu 185d71ae5a4SJacob 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) 186d71ae5a4SJacob Faibussowitsch { 1873d177a5cSEmil Constantinescu PetscFunctionBegin; 1883d177a5cSEmil Constantinescu PetscValidHeaderSpecific(adapt, TSGLLEADAPT_CLASSID, 1); 1894f572ea9SToby Isaac PetscAssertPointer(orders, 3); 1904f572ea9SToby Isaac PetscAssertPointer(errors, 4); 1914f572ea9SToby Isaac PetscAssertPointer(cost, 5); 1924f572ea9SToby Isaac PetscAssertPointer(next_sc, 9); 1934f572ea9SToby Isaac PetscAssertPointer(next_h, 10); 1944f572ea9SToby Isaac PetscAssertPointer(finish, 11); 195dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, n, orders, errors, cost, cur, h, tleft, next_sc, next_h, finish); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1973d177a5cSEmil Constantinescu } 1983d177a5cSEmil Constantinescu 199d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm, TSGLLEAdapt *inadapt) 200d71ae5a4SJacob Faibussowitsch { 2013d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 2023d177a5cSEmil Constantinescu 2033d177a5cSEmil Constantinescu PetscFunctionBegin; 2043d177a5cSEmil Constantinescu *inadapt = NULL; 2059566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSGLLEADAPT_CLASSID, "TSGLLEAdapt", "General Linear adaptivity", "TS", comm, TSGLLEAdaptDestroy, TSGLLEAdaptView)); 2063d177a5cSEmil Constantinescu *inadapt = adapt; 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2083d177a5cSEmil Constantinescu } 2093d177a5cSEmil Constantinescu 2103d177a5cSEmil Constantinescu /* 2110e3d61c9SBarry Smith Implementations 2123d177a5cSEmil Constantinescu */ 2133d177a5cSEmil Constantinescu 214d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) 215d71ae5a4SJacob Faibussowitsch { 2163d177a5cSEmil Constantinescu PetscFunctionBegin; 2179566063dSJacob Faibussowitsch PetscCall(PetscFree(adapt->data)); 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2193d177a5cSEmil Constantinescu } 2203d177a5cSEmil Constantinescu 2213d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */ 2223d177a5cSEmil Constantinescu typedef struct { 2233d177a5cSEmil Constantinescu PetscInt scheme; 2243d177a5cSEmil Constantinescu PetscReal h; 2253d177a5cSEmil Constantinescu } TSGLLEAdapt_None; 2263d177a5cSEmil Constantinescu 227d71ae5a4SJacob 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) 228d71ae5a4SJacob Faibussowitsch { 2293d177a5cSEmil Constantinescu PetscFunctionBegin; 2303d177a5cSEmil Constantinescu *next_sc = cur; 2313d177a5cSEmil Constantinescu *next_h = h; 2323d177a5cSEmil Constantinescu if (*next_h > tleft) { 2333d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2343d177a5cSEmil Constantinescu *next_h = tleft; 2353d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2373d177a5cSEmil Constantinescu } 2383d177a5cSEmil Constantinescu 239d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) 240d71ae5a4SJacob Faibussowitsch { 2413d177a5cSEmil Constantinescu TSGLLEAdapt_None *a; 2423d177a5cSEmil Constantinescu 2433d177a5cSEmil Constantinescu PetscFunctionBegin; 2444dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 2453d177a5cSEmil Constantinescu adapt->data = (void *)a; 2463d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_None; 2473d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2493d177a5cSEmil Constantinescu } 2503d177a5cSEmil Constantinescu 2513d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */ 2523d177a5cSEmil Constantinescu typedef struct { 2533d177a5cSEmil Constantinescu PetscReal desired_h; 2543d177a5cSEmil Constantinescu } TSGLLEAdapt_Size; 2553d177a5cSEmil Constantinescu 256d71ae5a4SJacob 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) 257d71ae5a4SJacob Faibussowitsch { 2583d177a5cSEmil Constantinescu TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size *)adapt->data; 2593d177a5cSEmil Constantinescu PetscReal dec = 0.2, inc = 5.0, safe = 0.9, optimal, last_desired_h; 2603d177a5cSEmil Constantinescu 2613d177a5cSEmil Constantinescu PetscFunctionBegin; 2623d177a5cSEmil Constantinescu *next_sc = cur; 2633d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[cur], (PetscReal)-1. / (safe * orders[cur])); 2643d177a5cSEmil Constantinescu /* Step sizes oscillate when there is no smoothing. Here we use a geometric mean of the current step size and the 2653d177a5cSEmil Constantinescu * one that would have been taken (without smoothing) on the last step. */ 2663d177a5cSEmil Constantinescu last_desired_h = sz->desired_h; 2673d177a5cSEmil Constantinescu sz->desired_h = h * PetscMax(dec, PetscMin(inc, optimal)); /* Trim to [dec,inc] */ 2683d177a5cSEmil Constantinescu 2693d177a5cSEmil Constantinescu /* Normally only happens on the first step */ 2703d177a5cSEmil Constantinescu if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h); 2713d177a5cSEmil Constantinescu else *next_h = sz->desired_h; 2723d177a5cSEmil Constantinescu 2733d177a5cSEmil Constantinescu if (*next_h > tleft) { 2743d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2753d177a5cSEmil Constantinescu *next_h = tleft; 2763d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2783d177a5cSEmil Constantinescu } 2793d177a5cSEmil Constantinescu 280d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) 281d71ae5a4SJacob Faibussowitsch { 2823d177a5cSEmil Constantinescu TSGLLEAdapt_Size *a; 2833d177a5cSEmil Constantinescu 2843d177a5cSEmil Constantinescu PetscFunctionBegin; 2854dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 2863d177a5cSEmil Constantinescu adapt->data = (void *)a; 2873d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Size; 2883d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2903d177a5cSEmil Constantinescu } 2913d177a5cSEmil Constantinescu 2923d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */ 2933d177a5cSEmil Constantinescu typedef struct { 2943d177a5cSEmil Constantinescu PetscInt count_at_order; 2953d177a5cSEmil Constantinescu PetscReal desired_h; 2963d177a5cSEmil Constantinescu } TSGLLEAdapt_Both; 2973d177a5cSEmil Constantinescu 298d71ae5a4SJacob 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) 299d71ae5a4SJacob Faibussowitsch { 3003d177a5cSEmil Constantinescu TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both *)adapt->data; 3013d177a5cSEmil Constantinescu PetscReal dec = 0.2, inc = 5.0, safe = 0.9; 3029371c9d4SSatish Balay struct { 3039371c9d4SSatish Balay PetscInt id; 3049371c9d4SSatish Balay PetscReal h, eff; 3059371c9d4SSatish Balay } best = {-1, 0, 0}, trial = {-1, 0, 0}, current = {-1, 0, 0}; 3063d177a5cSEmil Constantinescu PetscInt i; 3073d177a5cSEmil Constantinescu 3083d177a5cSEmil Constantinescu PetscFunctionBegin; 3093d177a5cSEmil Constantinescu for (i = 0; i < n; i++) { 3103d177a5cSEmil Constantinescu PetscReal optimal; 3113d177a5cSEmil Constantinescu trial.id = i; 3123d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[i], (PetscReal)-1. / (safe * orders[i])); 3133d177a5cSEmil Constantinescu trial.h = h * optimal; 3143d177a5cSEmil Constantinescu trial.eff = trial.h / cost[i]; 3159566063dSJacob Faibussowitsch if (trial.eff > best.eff) PetscCall(PetscArraycpy(&best, &trial, 1)); 3169566063dSJacob Faibussowitsch if (i == cur) PetscCall(PetscArraycpy(¤t, &trial, 1)); 3173d177a5cSEmil Constantinescu } 3183d177a5cSEmil Constantinescu /* Only switch orders if the scheme offers significant benefits over the current one. 3193d177a5cSEmil Constantinescu When the scheme is not changing, only change step size if it offers significant benefits. */ 3203d177a5cSEmil Constantinescu if (best.eff < 1.2 * current.eff || both->count_at_order < orders[cur] + 2) { 3213d177a5cSEmil Constantinescu PetscReal last_desired_h; 3223d177a5cSEmil Constantinescu *next_sc = current.id; 3233d177a5cSEmil Constantinescu last_desired_h = both->desired_h; 3243d177a5cSEmil Constantinescu both->desired_h = PetscMax(h * dec, PetscMin(h * inc, current.h)); 3259371c9d4SSatish Balay *next_h = (both->count_at_order > 0) ? PetscSqrtReal(last_desired_h * both->desired_h) : both->desired_h; 3263d177a5cSEmil Constantinescu both->count_at_order++; 3273d177a5cSEmil Constantinescu } else { 3283d177a5cSEmil Constantinescu PetscReal rat = cost[best.id] / cost[cur]; 3293d177a5cSEmil Constantinescu *next_sc = best.id; 3303d177a5cSEmil Constantinescu *next_h = PetscMax(h * rat * dec, PetscMin(h * rat * inc, best.h)); 3313d177a5cSEmil Constantinescu both->count_at_order = 0; 3323d177a5cSEmil Constantinescu both->desired_h = best.h; 3333d177a5cSEmil Constantinescu } 3343d177a5cSEmil Constantinescu 3353d177a5cSEmil Constantinescu if (*next_h > tleft) { 3363d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 3373d177a5cSEmil Constantinescu *next_h = tleft; 3383d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3403d177a5cSEmil Constantinescu } 3413d177a5cSEmil Constantinescu 342d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) 343d71ae5a4SJacob Faibussowitsch { 3443d177a5cSEmil Constantinescu TSGLLEAdapt_Both *a; 3453d177a5cSEmil Constantinescu 3463d177a5cSEmil Constantinescu PetscFunctionBegin; 3474dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 3483d177a5cSEmil Constantinescu adapt->data = (void *)a; 3493d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Both; 3503d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3523d177a5cSEmil Constantinescu } 353