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 26bcf0153eSBarry Smith TSGLLEAdaptRegister - adds a `TSGLLEAdapt` implementation 273d177a5cSEmil Constantinescu 283d177a5cSEmil Constantinescu Not Collective 293d177a5cSEmil Constantinescu 303d177a5cSEmil Constantinescu Input Parameters: 3120f4b53cSBarry Smith + sname - name of user-defined adaptivity scheme 3220f4b53cSBarry Smith - function - routine to create method context 3320f4b53cSBarry Smith 3420f4b53cSBarry Smith Level: advanced 353d177a5cSEmil Constantinescu 36bcf0153eSBarry Smith Note: 37bcf0153eSBarry Smith `TSGLLEAdaptRegister()` may be called multiple times to add several user-defined families. 383d177a5cSEmil Constantinescu 393d177a5cSEmil Constantinescu Sample usage: 403d177a5cSEmil Constantinescu .vb 413d177a5cSEmil Constantinescu TSGLLEAdaptRegister("my_scheme", MySchemeCreate); 423d177a5cSEmil Constantinescu .ve 433d177a5cSEmil Constantinescu 443d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 453d177a5cSEmil Constantinescu $ TSGLLEAdaptSetType(ts, "my_scheme") 463d177a5cSEmil Constantinescu or at runtime via the option 473d177a5cSEmil Constantinescu $ -ts_adapt_type my_scheme 483d177a5cSEmil Constantinescu 49*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLEAdapt`, `TSGLLEAdaptRegisterAll()` 503d177a5cSEmil Constantinescu @*/ 51d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegister(const char sname[], PetscErrorCode (*function)(TSGLLEAdapt)) 52d71ae5a4SJacob Faibussowitsch { 533d177a5cSEmil Constantinescu PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptInitializePackage()); 559566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEAdaptList, sname, function)); 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 573d177a5cSEmil Constantinescu } 583d177a5cSEmil Constantinescu 593d177a5cSEmil Constantinescu /*@C 60bcf0153eSBarry Smith TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in `TSGLLEAdapt` 613d177a5cSEmil Constantinescu 623d177a5cSEmil Constantinescu Not Collective 633d177a5cSEmil Constantinescu 643d177a5cSEmil Constantinescu Level: advanced 653d177a5cSEmil Constantinescu 66*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdapt`, `TSGLLE`, `TSGLLEAdaptRegisterDestroy()` 673d177a5cSEmil Constantinescu @*/ 68d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptRegisterAll(void) 69d71ae5a4SJacob Faibussowitsch { 703d177a5cSEmil Constantinescu PetscFunctionBegin; 713ba16761SJacob Faibussowitsch if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 723d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_TRUE; 739566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_NONE, TSGLLEAdaptCreate_None)); 749566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_SIZE, TSGLLEAdaptCreate_Size)); 759566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegister(TSGLLEADAPT_BOTH, TSGLLEAdaptCreate_Both)); 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 773d177a5cSEmil Constantinescu } 783d177a5cSEmil Constantinescu 793d177a5cSEmil Constantinescu /*@C 80bcf0153eSBarry Smith TSGLLEFinalizePackage - This function destroys everything in the `TSGLLE` package. It is 81bcf0153eSBarry Smith called from `PetscFinalize()`. 823d177a5cSEmil Constantinescu 833d177a5cSEmil Constantinescu Level: developer 843d177a5cSEmil Constantinescu 85*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSGLLEAdapt`, `TSGLLEAdaptInitializePackage()` 863d177a5cSEmil Constantinescu @*/ 87d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptFinalizePackage(void) 88d71ae5a4SJacob Faibussowitsch { 893d177a5cSEmil Constantinescu PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEAdaptList)); 913d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_FALSE; 923d177a5cSEmil Constantinescu TSGLLEAdaptRegisterAllCalled = PETSC_FALSE; 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 943d177a5cSEmil Constantinescu } 953d177a5cSEmil Constantinescu 963d177a5cSEmil Constantinescu /*@C 97bcf0153eSBarry Smith TSGLLEAdaptInitializePackage - This function initializes everything in the `TSGLLEAdapt` package. It is 98bcf0153eSBarry Smith called from `TSInitializePackage()`. 993d177a5cSEmil Constantinescu 1003d177a5cSEmil Constantinescu Level: developer 1013d177a5cSEmil Constantinescu 102*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSGLLEAdapt`, `TSGLLEAdaptFinalizePackage()` 1033d177a5cSEmil Constantinescu @*/ 104d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptInitializePackage(void) 105d71ae5a4SJacob Faibussowitsch { 1063d177a5cSEmil Constantinescu PetscFunctionBegin; 1073ba16761SJacob Faibussowitsch if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 1083d177a5cSEmil Constantinescu TSGLLEAdaptPackageInitialized = PETSC_TRUE; 1099566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSGLLEAdapt", &TSGLLEADAPT_CLASSID)); 1109566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptRegisterAll()); 1119566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSGLLEAdaptFinalizePackage)); 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1133d177a5cSEmil Constantinescu } 1143d177a5cSEmil Constantinescu 115d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt adapt, TSGLLEAdaptType type) 116d71ae5a4SJacob Faibussowitsch { 1175f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSGLLEAdapt); 1183d177a5cSEmil Constantinescu 1193d177a5cSEmil Constantinescu PetscFunctionBegin; 1209566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEAdaptList, type, &r)); 1213c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAdapt type \"%s\" given", type); 122dbbe0bcdSBarry Smith if (((PetscObject)adapt)->type_name) PetscUseTypeMethod(adapt, destroy); 1239566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 1249566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1263d177a5cSEmil Constantinescu } 1273d177a5cSEmil Constantinescu 128d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt, const char prefix[]) 129d71ae5a4SJacob Faibussowitsch { 1303d177a5cSEmil Constantinescu PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1333d177a5cSEmil Constantinescu } 1343d177a5cSEmil Constantinescu 135d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt adapt, PetscViewer viewer) 136d71ae5a4SJacob Faibussowitsch { 1373d177a5cSEmil Constantinescu PetscBool iascii; 1383d177a5cSEmil Constantinescu 1393d177a5cSEmil Constantinescu PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1413d177a5cSEmil Constantinescu if (iascii) { 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 1433d177a5cSEmil Constantinescu if (adapt->ops->view) { 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 145dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, view, viewer); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1473d177a5cSEmil Constantinescu } 1483d177a5cSEmil Constantinescu } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1503d177a5cSEmil Constantinescu } 1513d177a5cSEmil Constantinescu 152d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *adapt) 153d71ae5a4SJacob Faibussowitsch { 1543d177a5cSEmil Constantinescu PetscFunctionBegin; 1553ba16761SJacob Faibussowitsch if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS); 1563d177a5cSEmil Constantinescu PetscValidHeaderSpecific(*adapt, TSGLLEADAPT_CLASSID, 1); 1579371c9d4SSatish Balay if (--((PetscObject)(*adapt))->refct > 0) { 1589371c9d4SSatish Balay *adapt = NULL; 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1609371c9d4SSatish Balay } 161dbbe0bcdSBarry Smith PetscTryTypeMethod((*adapt), destroy); 1629566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1643d177a5cSEmil Constantinescu } 1653d177a5cSEmil Constantinescu 166d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt adapt, PetscOptionItems *PetscOptionsObject) 167d71ae5a4SJacob Faibussowitsch { 1683d177a5cSEmil Constantinescu char type[256] = TSGLLEADAPT_BOTH; 1693d177a5cSEmil Constantinescu PetscBool flg; 1703d177a5cSEmil Constantinescu 1713d177a5cSEmil Constantinescu PetscFunctionBegin; 1723d177a5cSEmil Constantinescu /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this 1733d177a5cSEmil Constantinescu * function can only be called from inside TSSetFromOptions_GLLE() */ 174d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TSGLLE Adaptivity options"); 1759371c9d4SSatish 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)); 1769566063dSJacob Faibussowitsch if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSGLLEAdaptSetType(adapt, type)); 177dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 178d0609cedSBarry Smith PetscOptionsHeadEnd(); 1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1803d177a5cSEmil Constantinescu } 1813d177a5cSEmil Constantinescu 182d71ae5a4SJacob 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) 183d71ae5a4SJacob Faibussowitsch { 1843d177a5cSEmil Constantinescu PetscFunctionBegin; 1853d177a5cSEmil Constantinescu PetscValidHeaderSpecific(adapt, TSGLLEADAPT_CLASSID, 1); 1863d177a5cSEmil Constantinescu PetscValidIntPointer(orders, 3); 187dadcf809SJacob Faibussowitsch PetscValidRealPointer(errors, 4); 188dadcf809SJacob Faibussowitsch PetscValidRealPointer(cost, 5); 1893d177a5cSEmil Constantinescu PetscValidIntPointer(next_sc, 9); 190dadcf809SJacob Faibussowitsch PetscValidRealPointer(next_h, 10); 191064a246eSJacob Faibussowitsch PetscValidBoolPointer(finish, 11); 192dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, n, orders, errors, cost, cur, h, tleft, next_sc, next_h, finish); 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1943d177a5cSEmil Constantinescu } 1953d177a5cSEmil Constantinescu 196d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate(MPI_Comm comm, TSGLLEAdapt *inadapt) 197d71ae5a4SJacob Faibussowitsch { 1983d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 1993d177a5cSEmil Constantinescu 2003d177a5cSEmil Constantinescu PetscFunctionBegin; 2013d177a5cSEmil Constantinescu *inadapt = NULL; 2029566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSGLLEADAPT_CLASSID, "TSGLLEAdapt", "General Linear adaptivity", "TS", comm, TSGLLEAdaptDestroy, TSGLLEAdaptView)); 2033d177a5cSEmil Constantinescu *inadapt = adapt; 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2053d177a5cSEmil Constantinescu } 2063d177a5cSEmil Constantinescu 2073d177a5cSEmil Constantinescu /* 2080e3d61c9SBarry Smith Implementations 2093d177a5cSEmil Constantinescu */ 2103d177a5cSEmil Constantinescu 211d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt) 212d71ae5a4SJacob Faibussowitsch { 2133d177a5cSEmil Constantinescu PetscFunctionBegin; 2149566063dSJacob Faibussowitsch PetscCall(PetscFree(adapt->data)); 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2163d177a5cSEmil Constantinescu } 2173d177a5cSEmil Constantinescu 2183d177a5cSEmil Constantinescu /* -------------------------------- None ----------------------------------- */ 2193d177a5cSEmil Constantinescu typedef struct { 2203d177a5cSEmil Constantinescu PetscInt scheme; 2213d177a5cSEmil Constantinescu PetscReal h; 2223d177a5cSEmil Constantinescu } TSGLLEAdapt_None; 2233d177a5cSEmil Constantinescu 224d71ae5a4SJacob 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) 225d71ae5a4SJacob Faibussowitsch { 2263d177a5cSEmil Constantinescu PetscFunctionBegin; 2273d177a5cSEmil Constantinescu *next_sc = cur; 2283d177a5cSEmil Constantinescu *next_h = h; 2293d177a5cSEmil Constantinescu if (*next_h > tleft) { 2303d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2313d177a5cSEmil Constantinescu *next_h = tleft; 2323d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2343d177a5cSEmil Constantinescu } 2353d177a5cSEmil Constantinescu 236d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt adapt) 237d71ae5a4SJacob Faibussowitsch { 2383d177a5cSEmil Constantinescu TSGLLEAdapt_None *a; 2393d177a5cSEmil Constantinescu 2403d177a5cSEmil Constantinescu PetscFunctionBegin; 2414dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 2423d177a5cSEmil Constantinescu adapt->data = (void *)a; 2433d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_None; 2443d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2463d177a5cSEmil Constantinescu } 2473d177a5cSEmil Constantinescu 2483d177a5cSEmil Constantinescu /* -------------------------------- Size ----------------------------------- */ 2493d177a5cSEmil Constantinescu typedef struct { 2503d177a5cSEmil Constantinescu PetscReal desired_h; 2513d177a5cSEmil Constantinescu } TSGLLEAdapt_Size; 2523d177a5cSEmil Constantinescu 253d71ae5a4SJacob 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) 254d71ae5a4SJacob Faibussowitsch { 2553d177a5cSEmil Constantinescu TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size *)adapt->data; 2563d177a5cSEmil Constantinescu PetscReal dec = 0.2, inc = 5.0, safe = 0.9, optimal, last_desired_h; 2573d177a5cSEmil Constantinescu 2583d177a5cSEmil Constantinescu PetscFunctionBegin; 2593d177a5cSEmil Constantinescu *next_sc = cur; 2603d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[cur], (PetscReal)-1. / (safe * orders[cur])); 2613d177a5cSEmil Constantinescu /* Step sizes oscillate when there is no smoothing. Here we use a geometric mean of the current step size and the 2623d177a5cSEmil Constantinescu * one that would have been taken (without smoothing) on the last step. */ 2633d177a5cSEmil Constantinescu last_desired_h = sz->desired_h; 2643d177a5cSEmil Constantinescu sz->desired_h = h * PetscMax(dec, PetscMin(inc, optimal)); /* Trim to [dec,inc] */ 2653d177a5cSEmil Constantinescu 2663d177a5cSEmil Constantinescu /* Normally only happens on the first step */ 2673d177a5cSEmil Constantinescu if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h); 2683d177a5cSEmil Constantinescu else *next_h = sz->desired_h; 2693d177a5cSEmil Constantinescu 2703d177a5cSEmil Constantinescu if (*next_h > tleft) { 2713d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 2723d177a5cSEmil Constantinescu *next_h = tleft; 2733d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2753d177a5cSEmil Constantinescu } 2763d177a5cSEmil Constantinescu 277d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt) 278d71ae5a4SJacob Faibussowitsch { 2793d177a5cSEmil Constantinescu TSGLLEAdapt_Size *a; 2803d177a5cSEmil Constantinescu 2813d177a5cSEmil Constantinescu PetscFunctionBegin; 2824dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 2833d177a5cSEmil Constantinescu adapt->data = (void *)a; 2843d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Size; 2853d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2873d177a5cSEmil Constantinescu } 2883d177a5cSEmil Constantinescu 2893d177a5cSEmil Constantinescu /* -------------------------------- Both ----------------------------------- */ 2903d177a5cSEmil Constantinescu typedef struct { 2913d177a5cSEmil Constantinescu PetscInt count_at_order; 2923d177a5cSEmil Constantinescu PetscReal desired_h; 2933d177a5cSEmil Constantinescu } TSGLLEAdapt_Both; 2943d177a5cSEmil Constantinescu 295d71ae5a4SJacob 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) 296d71ae5a4SJacob Faibussowitsch { 2973d177a5cSEmil Constantinescu TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both *)adapt->data; 2983d177a5cSEmil Constantinescu PetscReal dec = 0.2, inc = 5.0, safe = 0.9; 2999371c9d4SSatish Balay struct { 3009371c9d4SSatish Balay PetscInt id; 3019371c9d4SSatish Balay PetscReal h, eff; 3029371c9d4SSatish Balay } best = {-1, 0, 0}, trial = {-1, 0, 0}, current = {-1, 0, 0}; 3033d177a5cSEmil Constantinescu PetscInt i; 3043d177a5cSEmil Constantinescu 3053d177a5cSEmil Constantinescu PetscFunctionBegin; 3063d177a5cSEmil Constantinescu for (i = 0; i < n; i++) { 3073d177a5cSEmil Constantinescu PetscReal optimal; 3083d177a5cSEmil Constantinescu trial.id = i; 3093d177a5cSEmil Constantinescu optimal = PetscPowReal((PetscReal)errors[i], (PetscReal)-1. / (safe * orders[i])); 3103d177a5cSEmil Constantinescu trial.h = h * optimal; 3113d177a5cSEmil Constantinescu trial.eff = trial.h / cost[i]; 3129566063dSJacob Faibussowitsch if (trial.eff > best.eff) PetscCall(PetscArraycpy(&best, &trial, 1)); 3139566063dSJacob Faibussowitsch if (i == cur) PetscCall(PetscArraycpy(¤t, &trial, 1)); 3143d177a5cSEmil Constantinescu } 3153d177a5cSEmil Constantinescu /* Only switch orders if the scheme offers significant benefits over the current one. 3163d177a5cSEmil Constantinescu When the scheme is not changing, only change step size if it offers significant benefits. */ 3173d177a5cSEmil Constantinescu if (best.eff < 1.2 * current.eff || both->count_at_order < orders[cur] + 2) { 3183d177a5cSEmil Constantinescu PetscReal last_desired_h; 3193d177a5cSEmil Constantinescu *next_sc = current.id; 3203d177a5cSEmil Constantinescu last_desired_h = both->desired_h; 3213d177a5cSEmil Constantinescu both->desired_h = PetscMax(h * dec, PetscMin(h * inc, current.h)); 3229371c9d4SSatish Balay *next_h = (both->count_at_order > 0) ? PetscSqrtReal(last_desired_h * both->desired_h) : both->desired_h; 3233d177a5cSEmil Constantinescu both->count_at_order++; 3243d177a5cSEmil Constantinescu } else { 3253d177a5cSEmil Constantinescu PetscReal rat = cost[best.id] / cost[cur]; 3263d177a5cSEmil Constantinescu *next_sc = best.id; 3273d177a5cSEmil Constantinescu *next_h = PetscMax(h * rat * dec, PetscMin(h * rat * inc, best.h)); 3283d177a5cSEmil Constantinescu both->count_at_order = 0; 3293d177a5cSEmil Constantinescu both->desired_h = best.h; 3303d177a5cSEmil Constantinescu } 3313d177a5cSEmil Constantinescu 3323d177a5cSEmil Constantinescu if (*next_h > tleft) { 3333d177a5cSEmil Constantinescu *finish = PETSC_TRUE; 3343d177a5cSEmil Constantinescu *next_h = tleft; 3353d177a5cSEmil Constantinescu } else *finish = PETSC_FALSE; 3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3373d177a5cSEmil Constantinescu } 3383d177a5cSEmil Constantinescu 339d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt) 340d71ae5a4SJacob Faibussowitsch { 3413d177a5cSEmil Constantinescu TSGLLEAdapt_Both *a; 3423d177a5cSEmil Constantinescu 3433d177a5cSEmil Constantinescu PetscFunctionBegin; 3444dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 3453d177a5cSEmil Constantinescu adapt->data = (void *)a; 3463d177a5cSEmil Constantinescu adapt->ops->choose = TSGLLEAdaptChoose_Both; 3473d177a5cSEmil Constantinescu adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree; 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3493d177a5cSEmil Constantinescu } 350