1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 284df9cb4SJed Brown 31b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID; 41b9b13dfSLisandro Dalcin 5140e18c1SBarry Smith static PetscFunctionList TSAdaptList; 684df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 784df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 884df9cb4SJed Brown 98cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 101566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 11cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt); 128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 131917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 14d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt); 1584df9cb4SJed Brown 1684df9cb4SJed Brown /*@C 171c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation 181c84c290SBarry Smith 19cc4c1da9SBarry Smith Not Collective, No Fortran Support 201c84c290SBarry Smith 211c84c290SBarry Smith Input Parameters: 222fe279fdSBarry Smith + sname - name of user-defined adaptivity scheme 232fe279fdSBarry Smith - function - routine to create method context 241c84c290SBarry Smith 25bcf0153eSBarry Smith Level: advanced 26bcf0153eSBarry Smith 271c84c290SBarry Smith Notes: 28bcf0153eSBarry Smith `TSAdaptRegister()` may be called multiple times to add several user-defined families. 291c84c290SBarry Smith 30b43aa488SJacob Faibussowitsch Example Usage: 311c84c290SBarry Smith .vb 32bdf89e91SBarry Smith TSAdaptRegister("my_scheme", MySchemeCreate); 331c84c290SBarry Smith .ve 341c84c290SBarry Smith 351c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via 361c84c290SBarry Smith $ TSAdaptSetType(ts, "my_scheme") 371c84c290SBarry Smith or at runtime via the option 381c84c290SBarry Smith $ -ts_adapt_type my_scheme 3984df9cb4SJed Brown 401cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterAll()` 4184df9cb4SJed Brown @*/ 42d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt)) 43d71ae5a4SJacob Faibussowitsch { 4484df9cb4SJed Brown PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 469566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function)); 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4884df9cb4SJed Brown } 4984df9cb4SJed Brown 5084df9cb4SJed Brown /*@C 512fe279fdSBarry Smith TSAdaptRegisterAll - Registers all of the adaptivity schemes in `TSAdapt` 5284df9cb4SJed Brown 5384df9cb4SJed Brown Not Collective 5484df9cb4SJed Brown 5584df9cb4SJed Brown Level: advanced 5684df9cb4SJed Brown 571cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterDestroy()` 5884df9cb4SJed Brown @*/ 59d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegisterAll(void) 60d71ae5a4SJacob Faibussowitsch { 6184df9cb4SJed Brown PetscFunctionBegin; 623ba16761SJacob Faibussowitsch if (TSAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 630f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 649566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None)); 659566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic)); 669566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP)); 679566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL)); 689566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE)); 699566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History)); 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7184df9cb4SJed Brown } 7284df9cb4SJed Brown 7384df9cb4SJed Brown /*@C 742fe279fdSBarry Smith TSAdaptFinalizePackage - This function destroys everything in the `TS` package. It is 752fe279fdSBarry Smith called from `PetscFinalize()`. 7684df9cb4SJed Brown 7784df9cb4SJed Brown Level: developer 7884df9cb4SJed Brown 791cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()` 8084df9cb4SJed Brown @*/ 81d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptFinalizePackage(void) 82d71ae5a4SJacob Faibussowitsch { 8384df9cb4SJed Brown PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSAdaptList)); 8584df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 8684df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8884df9cb4SJed Brown } 8984df9cb4SJed Brown 9084df9cb4SJed Brown /*@C 912fe279fdSBarry Smith TSAdaptInitializePackage - This function initializes everything in the `TSAdapt` package. It is 922fe279fdSBarry Smith called from `TSInitializePackage()`. 9384df9cb4SJed Brown 9484df9cb4SJed Brown Level: developer 9584df9cb4SJed Brown 961cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()` 9784df9cb4SJed Brown @*/ 98d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptInitializePackage(void) 99d71ae5a4SJacob Faibussowitsch { 10084df9cb4SJed Brown PetscFunctionBegin; 1013ba16761SJacob Faibussowitsch if (TSAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 10284df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 1039566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID)); 1049566063dSJacob Faibussowitsch PetscCall(TSAdaptRegisterAll()); 1059566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage)); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10784df9cb4SJed Brown } 10884df9cb4SJed Brown 109cc4c1da9SBarry Smith /*@ 110bcf0153eSBarry Smith TSAdaptSetType - sets the approach used for the error adapter 1117eef6055SBarry Smith 11220f4b53cSBarry Smith Logicially Collective 1137eef6055SBarry Smith 114d8d19677SJose E. Roman Input Parameters: 11558b119d1SBarry Smith + adapt - the `TS` adapter, most likely obtained with `TSGetAdapt()` 116bcf0153eSBarry Smith - type - one of the `TSAdaptType` 1177eef6055SBarry Smith 118bcf0153eSBarry Smith Options Database Key: 119ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1207eef6055SBarry Smith 1217eef6055SBarry Smith Level: intermediate 1227eef6055SBarry Smith 123b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()` 1247eef6055SBarry Smith @*/ 125d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type) 126d71ae5a4SJacob Faibussowitsch { 127ef749922SLisandro Dalcin PetscBool match; 1285f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSAdapt); 12984df9cb4SJed Brown 13084df9cb4SJed Brown PetscFunctionBegin; 1314782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1324f572ea9SToby Isaac PetscAssertPointer(type, 2); 1339566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match)); 1343ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1359566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSAdaptList, type, &r)); 1366adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type); 137dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, destroy); 1389566063dSJacob Faibussowitsch PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps))); 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 1409566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14284df9cb4SJed Brown } 14384df9cb4SJed Brown 144cc4c1da9SBarry Smith /*@ 145a1cb98faSBarry Smith TSAdaptGetType - gets the `TS` adapter method type (as a string). 146d0288e4fSLisandro Dalcin 147d0288e4fSLisandro Dalcin Not Collective 148d0288e4fSLisandro Dalcin 149d0288e4fSLisandro Dalcin Input Parameter: 150a1cb98faSBarry Smith . adapt - The `TS` adapter, most likely obtained with `TSGetAdapt()` 151d0288e4fSLisandro Dalcin 152d0288e4fSLisandro Dalcin Output Parameter: 153a1cb98faSBarry Smith . type - The name of `TS` adapter method 154d0288e4fSLisandro Dalcin 155d0288e4fSLisandro Dalcin Level: intermediate 156d0288e4fSLisandro Dalcin 157a1cb98faSBarry Smith .seealso: `TSAdapt`, `TSAdaptType`, `TSAdaptSetType()` 158d0288e4fSLisandro Dalcin @*/ 159d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type) 160d71ae5a4SJacob Faibussowitsch { 161d0288e4fSLisandro Dalcin PetscFunctionBegin; 162d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1634f572ea9SToby Isaac PetscAssertPointer(type, 2); 164d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166d0288e4fSLisandro Dalcin } 167d0288e4fSLisandro Dalcin 168d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[]) 169d71ae5a4SJacob Faibussowitsch { 17084df9cb4SJed Brown PetscFunctionBegin; 1714782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17484df9cb4SJed Brown } 17584df9cb4SJed Brown 176ffeef943SBarry Smith /*@ 1772fe279fdSBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with `TSAdaptView()`. 178ad6bc421SBarry Smith 179c3339decSBarry Smith Collective 180ad6bc421SBarry Smith 181ad6bc421SBarry Smith Input Parameters: 182b43aa488SJacob Faibussowitsch + adapt - the newly loaded `TSAdapt`, this needs to have been created with `TSAdaptCreate()` or 183bcf0153eSBarry Smith some related function before a call to `TSAdaptLoad()`. 184bcf0153eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or 185bcf0153eSBarry Smith HDF5 file viewer, obtained from `PetscViewerHDF5Open()` 186ad6bc421SBarry Smith 187ad6bc421SBarry Smith Level: intermediate 188ad6bc421SBarry Smith 189bcf0153eSBarry Smith Note: 190bcf0153eSBarry Smith The type is determined by the data in the file, any type set into the `TSAdapt` before this call is ignored. 191ad6bc421SBarry Smith 1921cc06b55SBarry Smith .seealso: [](ch_ts), `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()`, `TSAdapt` 193ad6bc421SBarry Smith @*/ 194d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer) 195d71ae5a4SJacob Faibussowitsch { 196ad6bc421SBarry Smith PetscBool isbinary; 197ad6bc421SBarry Smith char type[256]; 198ad6bc421SBarry Smith 199ad6bc421SBarry Smith PetscFunctionBegin; 2004782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 201ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 2033c633725SBarry Smith PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 204ad6bc421SBarry Smith 2059566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 2069566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, type)); 207dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, load, viewer); 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209ad6bc421SBarry Smith } 210ad6bc421SBarry Smith 211d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer) 212d71ae5a4SJacob Faibussowitsch { 2131c167fc2SEmil Constantinescu PetscBool iascii, isbinary, isnone, isglee; 21484df9cb4SJed Brown 21584df9cb4SJed Brown PetscFunctionBegin; 2164782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 2179566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer)); 2184782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2194782b174SLisandro Dalcin PetscCheckSameComm(adapt, 1, viewer, 2); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 22284df9cb4SJed Brown if (iascii) { 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee)); 2261917a363SLisandro Dalcin if (!isnone) { 2279566063dSJacob Faibussowitsch if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, " always accepting steps\n")); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " safety factor %g\n", (double)adapt->safety)); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " extra safety factor after step rejection %g\n", (double)adapt->reject_safety)); 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest increase %g\n", (double)adapt->clip[1])); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest decrease %g\n", (double)adapt->clip[0])); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum allowed timestep %g\n", (double)adapt->dt_max)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " minimum allowed timestep %g\n", (double)adapt->dt_min)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max)); 2351c167fc2SEmil Constantinescu } 2361c167fc2SEmil Constantinescu if (isglee) { 2371c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses local error control\n")); 2391c167fc2SEmil Constantinescu } else { 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses global error control\n")); 2411c167fc2SEmil Constantinescu } 2421917a363SLisandro Dalcin } 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 244dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, view, viewer); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 246ad6bc421SBarry Smith } else if (isbinary) { 247ad6bc421SBarry Smith char type[256]; 248ad6bc421SBarry Smith 249ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 2509566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256)); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 252dbbe0bcdSBarry Smith } else PetscTryTypeMethod(adapt, view, viewer); 2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25484df9cb4SJed Brown } 25584df9cb4SJed Brown 256099af0a3SLisandro Dalcin /*@ 257bcf0153eSBarry Smith TSAdaptReset - Resets a `TSAdapt` context to its defaults 258099af0a3SLisandro Dalcin 259c3339decSBarry Smith Collective 260099af0a3SLisandro Dalcin 261099af0a3SLisandro Dalcin Input Parameter: 262bcf0153eSBarry Smith . adapt - the `TSAdapt` context obtained from `TSGetAdapt()` or `TSAdaptCreate()` 263099af0a3SLisandro Dalcin 264099af0a3SLisandro Dalcin Level: developer 265099af0a3SLisandro Dalcin 2661cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdapt`, `TSAdaptCreate()`, `TSAdaptDestroy()` 267099af0a3SLisandro Dalcin @*/ 268d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptReset(TSAdapt adapt) 269d71ae5a4SJacob Faibussowitsch { 270099af0a3SLisandro Dalcin PetscFunctionBegin; 271099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 272dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, reset); 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274099af0a3SLisandro Dalcin } 275099af0a3SLisandro Dalcin 276d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 277d71ae5a4SJacob Faibussowitsch { 27884df9cb4SJed Brown PetscFunctionBegin; 2793ba16761SJacob Faibussowitsch if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS); 28084df9cb4SJed Brown PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1); 281f4f49eeaSPierre Jolivet if (--((PetscObject)*adapt)->refct > 0) { 2829371c9d4SSatish Balay *adapt = NULL; 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2849371c9d4SSatish Balay } 285099af0a3SLisandro Dalcin 2869566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt)); 287099af0a3SLisandro Dalcin 288f4f49eeaSPierre Jolivet PetscTryTypeMethod(*adapt, destroy); 2899566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 2909566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29284df9cb4SJed Brown } 29384df9cb4SJed Brown 2941c3436cfSJed Brown /*@ 2951c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2961c3436cfSJed Brown 297c3339decSBarry Smith Collective 2981c3436cfSJed Brown 2994165533cSJose E. Roman Input Parameters: 3001c3436cfSJed Brown + adapt - adaptive controller context 301bcf0153eSBarry Smith - flg - `PETSC_TRUE` to active a monitor, `PETSC_FALSE` to disable 3021c3436cfSJed Brown 303bcf0153eSBarry Smith Options Database Key: 304ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 305bf997491SLisandro Dalcin 3061c3436cfSJed Brown Level: intermediate 3071c3436cfSJed Brown 3081cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 3091c3436cfSJed Brown @*/ 310d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg) 311d71ae5a4SJacob Faibussowitsch { 3121c3436cfSJed Brown PetscFunctionBegin; 3134782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3144782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flg, 2); 3151c3436cfSJed Brown if (flg) { 3169566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor)); 3171c3436cfSJed Brown } else { 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor)); 3191c3436cfSJed Brown } 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3211c3436cfSJed Brown } 3221c3436cfSJed Brown 3230873808bSJed Brown /*@C 324bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3250873808bSJed Brown 3262fe279fdSBarry Smith Logically Collective 3270873808bSJed Brown 3284165533cSJose E. Roman Input Parameters: 3290873808bSJed Brown + adapt - adaptive controller context 3300873808bSJed Brown - func - stage check function 3310873808bSJed Brown 332b43aa488SJacob Faibussowitsch Calling sequence: 3330873808bSJed Brown + adapt - adaptive controller context 3340873808bSJed Brown . ts - time stepping context 33514d0ab18SJacob Faibussowitsch . t - current time 33614d0ab18SJacob Faibussowitsch . Y - current solution vector 3370873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3380873808bSJed Brown 3390873808bSJed Brown Level: advanced 3400873808bSJed Brown 3411cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 3420873808bSJed Brown @*/ 34314d0ab18SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)) 344d71ae5a4SJacob Faibussowitsch { 3450873808bSJed Brown PetscFunctionBegin; 3460873808bSJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 34768ae941aSLisandro Dalcin adapt->checkstage = func; 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3490873808bSJed Brown } 3500873808bSJed Brown 3511c3436cfSJed Brown /*@ 352bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 353bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 354bf997491SLisandro Dalcin 3552fe279fdSBarry Smith Logically Collective 356bf997491SLisandro Dalcin 3574165533cSJose E. Roman Input Parameters: 358bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 359bf997491SLisandro Dalcin - flag - whether to always accept steps 360bf997491SLisandro Dalcin 361bcf0153eSBarry Smith Options Database Key: 362ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 363bf997491SLisandro Dalcin 364bf997491SLisandro Dalcin Level: intermediate 365bf997491SLisandro Dalcin 3661cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 367bf997491SLisandro Dalcin @*/ 368d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag) 369d71ae5a4SJacob Faibussowitsch { 370bf997491SLisandro Dalcin PetscFunctionBegin; 371bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 372bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flag, 2); 373bf997491SLisandro Dalcin adapt->always_accept = flag; 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 375bf997491SLisandro Dalcin } 376bf997491SLisandro Dalcin 377bf997491SLisandro Dalcin /*@ 378bcf0153eSBarry Smith TSAdaptSetSafety - Set safety factors for time step adaptor 3791c3436cfSJed Brown 3802fe279fdSBarry Smith Logically Collective 3811917a363SLisandro Dalcin 3824165533cSJose E. Roman Input Parameters: 3831917a363SLisandro Dalcin + adapt - adaptive controller context 3841917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 385ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 3861917a363SLisandro Dalcin 3871917a363SLisandro Dalcin Options Database Keys: 388ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 389ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 3901917a363SLisandro Dalcin 3911917a363SLisandro Dalcin Level: intermediate 3921917a363SLisandro Dalcin 39309cb0f53SBarry Smith Note: 39409cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 39509cb0f53SBarry Smith 39609cb0f53SBarry Smith Fortran Note: 39709cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 39809cb0f53SBarry Smith 3991cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 4001917a363SLisandro Dalcin @*/ 401d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety) 402d71ae5a4SJacob Faibussowitsch { 4031917a363SLisandro Dalcin PetscFunctionBegin; 4041917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4051917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, safety, 2); 4061917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, reject_safety, 3); 40709cb0f53SBarry Smith PetscCheck(safety == (PetscReal)PETSC_CURRENT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety); 40809cb0f53SBarry Smith PetscCheck(safety == (PetscReal)PETSC_CURRENT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety); 40909cb0f53SBarry Smith PetscCheck(reject_safety == (PetscReal)PETSC_CURRENT || reject_safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be non negative", (double)reject_safety); 41009cb0f53SBarry Smith PetscCheck(reject_safety == (PetscReal)PETSC_CURRENT || reject_safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be less than one", (double)reject_safety); 41109cb0f53SBarry Smith if (safety != (PetscReal)PETSC_CURRENT) adapt->safety = safety; 41209cb0f53SBarry Smith if (reject_safety != (PetscReal)PETSC_CURRENT) adapt->reject_safety = reject_safety; 4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4141917a363SLisandro Dalcin } 4151917a363SLisandro Dalcin 4161917a363SLisandro Dalcin /*@ 417bcf0153eSBarry Smith TSAdaptGetSafety - Get safety factors for time step adapter 4181917a363SLisandro Dalcin 4191917a363SLisandro Dalcin Not Collective 4201917a363SLisandro Dalcin 4214165533cSJose E. Roman Input Parameter: 4221917a363SLisandro Dalcin . adapt - adaptive controller context 4231917a363SLisandro Dalcin 4244165533cSJose E. Roman Output Parameters: 425b43aa488SJacob Faibussowitsch + safety - safety factor relative to target error/stability goal 426b43aa488SJacob Faibussowitsch - reject_safety - extra safety factor to apply if the last step was rejected 4271917a363SLisandro Dalcin 4281917a363SLisandro Dalcin Level: intermediate 4291917a363SLisandro Dalcin 4301cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 4311917a363SLisandro Dalcin @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety) 433d71ae5a4SJacob Faibussowitsch { 4341917a363SLisandro Dalcin PetscFunctionBegin; 4351917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4364f572ea9SToby Isaac if (safety) PetscAssertPointer(safety, 2); 4374f572ea9SToby Isaac if (reject_safety) PetscAssertPointer(reject_safety, 3); 4381917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4391917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4411917a363SLisandro Dalcin } 4421917a363SLisandro Dalcin 4431917a363SLisandro Dalcin /*@ 444bcf0153eSBarry Smith TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms 445bcf0153eSBarry Smith for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components. 44676cddca1SEmil Constantinescu 4472fe279fdSBarry Smith Logically Collective 44876cddca1SEmil Constantinescu 4494165533cSJose E. Roman Input Parameters: 45076cddca1SEmil Constantinescu + adapt - adaptive controller context 45176cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 45276cddca1SEmil Constantinescu 453bcf0153eSBarry Smith Options Database Key: 45476cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 45576cddca1SEmil Constantinescu 45676cddca1SEmil Constantinescu Level: intermediate 45776cddca1SEmil Constantinescu 4581cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 45976cddca1SEmil Constantinescu @*/ 460d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore) 461d71ae5a4SJacob Faibussowitsch { 46276cddca1SEmil Constantinescu PetscFunctionBegin; 46376cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 46476cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt, max_ignore, 2); 46576cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 4663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46776cddca1SEmil Constantinescu } 46876cddca1SEmil Constantinescu 46976cddca1SEmil Constantinescu /*@ 470bcf0153eSBarry Smith TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms 471bcf0153eSBarry Smith for time step adaptivity (in absolute value). 47276cddca1SEmil Constantinescu 47376cddca1SEmil Constantinescu Not Collective 47476cddca1SEmil Constantinescu 4754165533cSJose E. Roman Input Parameter: 47676cddca1SEmil Constantinescu . adapt - adaptive controller context 47776cddca1SEmil Constantinescu 4784165533cSJose E. Roman Output Parameter: 47976cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 48076cddca1SEmil Constantinescu 48176cddca1SEmil Constantinescu Level: intermediate 48276cddca1SEmil Constantinescu 4831cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 48476cddca1SEmil Constantinescu @*/ 485d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore) 486d71ae5a4SJacob Faibussowitsch { 48776cddca1SEmil Constantinescu PetscFunctionBegin; 48876cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4894f572ea9SToby Isaac PetscAssertPointer(max_ignore, 2); 49076cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49276cddca1SEmil Constantinescu } 49376cddca1SEmil Constantinescu 49476cddca1SEmil Constantinescu /*@ 495bcf0153eSBarry Smith TSAdaptSetClip - Sets the admissible decrease/increase factor in step size in the time step adapter 4961917a363SLisandro Dalcin 497c3339decSBarry Smith Logically collective 4981917a363SLisandro Dalcin 4994165533cSJose E. Roman Input Parameters: 5001917a363SLisandro Dalcin + adapt - adaptive controller context 5011917a363SLisandro Dalcin . low - admissible decrease factor 502ec18b7bcSLisandro Dalcin - high - admissible increase factor 5031917a363SLisandro Dalcin 504bcf0153eSBarry Smith Options Database Key: 505ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5061917a363SLisandro Dalcin 5071917a363SLisandro Dalcin Level: intermediate 5081917a363SLisandro Dalcin 50909cb0f53SBarry Smith Note: 51009cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 51109cb0f53SBarry Smith 51209cb0f53SBarry Smith Fortran Note: 51309cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 51409cb0f53SBarry Smith 5151cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 5161917a363SLisandro Dalcin @*/ 517d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high) 518d71ae5a4SJacob Faibussowitsch { 5191917a363SLisandro Dalcin PetscFunctionBegin; 5201917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 5211917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, low, 2); 5221917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, high, 3); 52309cb0f53SBarry Smith PetscCheck(low == (PetscReal)PETSC_CURRENT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low); 52409cb0f53SBarry Smith PetscCheck(low == (PetscReal)PETSC_CURRENT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low); 52509cb0f53SBarry Smith PetscCheck(high == (PetscReal)PETSC_CURRENT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high); 52609cb0f53SBarry Smith if (low != (PetscReal)PETSC_CURRENT) adapt->clip[0] = low; 52709cb0f53SBarry Smith if (high != (PetscReal)PETSC_CURRENT) adapt->clip[1] = high; 5283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5291917a363SLisandro Dalcin } 5301917a363SLisandro Dalcin 5311917a363SLisandro Dalcin /*@ 532bcf0153eSBarry Smith TSAdaptGetClip - Gets the admissible decrease/increase factor in step size in the time step adapter 5331917a363SLisandro Dalcin 5341917a363SLisandro Dalcin Not Collective 5351917a363SLisandro Dalcin 5364165533cSJose E. Roman Input Parameter: 5371917a363SLisandro Dalcin . adapt - adaptive controller context 5381917a363SLisandro Dalcin 5394165533cSJose E. Roman Output Parameters: 5401917a363SLisandro Dalcin + low - optional, admissible decrease factor 5411917a363SLisandro Dalcin - high - optional, admissible increase factor 5421917a363SLisandro Dalcin 5431917a363SLisandro Dalcin Level: intermediate 5441917a363SLisandro Dalcin 5451cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 5461917a363SLisandro Dalcin @*/ 547d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high) 548d71ae5a4SJacob Faibussowitsch { 5491917a363SLisandro Dalcin PetscFunctionBegin; 5501917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 5514f572ea9SToby Isaac if (low) PetscAssertPointer(low, 2); 5524f572ea9SToby Isaac if (high) PetscAssertPointer(high, 3); 5531917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5541917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5561917a363SLisandro Dalcin } 5571917a363SLisandro Dalcin 5581917a363SLisandro Dalcin /*@ 559bcf0153eSBarry Smith TSAdaptSetScaleSolveFailed - Scale step size by this factor if solve fails 56062c23b28SMark 56120f4b53cSBarry Smith Logically Collective 56262c23b28SMark 5634165533cSJose E. Roman Input Parameters: 56462c23b28SMark + adapt - adaptive controller context 565ee300463SSatish Balay - scale - scale 56662c23b28SMark 567bcf0153eSBarry Smith Options Database Key: 56862c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 56962c23b28SMark 57062c23b28SMark Level: intermediate 57162c23b28SMark 5721cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 57362c23b28SMark @*/ 574d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale) 575d71ae5a4SJacob Faibussowitsch { 57662c23b28SMark PetscFunctionBegin; 57762c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 57862c23b28SMark PetscValidLogicalCollectiveReal(adapt, scale, 2); 57909cb0f53SBarry Smith PetscCheck(scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale); 58009cb0f53SBarry Smith PetscCheck(scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale); 58109cb0f53SBarry Smith adapt->scale_solve_failed = scale; 5823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58362c23b28SMark } 58462c23b28SMark 58562c23b28SMark /*@ 58662c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 58762c23b28SMark 58862c23b28SMark Not Collective 58962c23b28SMark 5904165533cSJose E. Roman Input Parameter: 59162c23b28SMark . adapt - adaptive controller context 59262c23b28SMark 5934165533cSJose E. Roman Output Parameter: 594ee300463SSatish Balay . scale - scale factor 59562c23b28SMark 59662c23b28SMark Level: intermediate 59762c23b28SMark 5981cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 59962c23b28SMark @*/ 600d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale) 601d71ae5a4SJacob Faibussowitsch { 60262c23b28SMark PetscFunctionBegin; 60362c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 6044f572ea9SToby Isaac if (scale) PetscAssertPointer(scale, 2); 60562c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 6063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60762c23b28SMark } 60862c23b28SMark 60962c23b28SMark /*@ 610bcf0153eSBarry Smith TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the time step controller 6111917a363SLisandro Dalcin 61220f4b53cSBarry Smith Logically Collective 6131c3436cfSJed Brown 6144165533cSJose E. Roman Input Parameters: 615bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 6161c3436cfSJed Brown . hmin - minimum time step 6171c3436cfSJed Brown - hmax - maximum time step 6181c3436cfSJed Brown 6191c3436cfSJed Brown Options Database Keys: 620ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 621ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 6221c3436cfSJed Brown 6231c3436cfSJed Brown Level: intermediate 6241c3436cfSJed Brown 62509cb0f53SBarry Smith Note: 62609cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 62709cb0f53SBarry Smith 62809cb0f53SBarry Smith Fortran Note: 62909cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 63009cb0f53SBarry Smith 6311cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 6321c3436cfSJed Brown @*/ 633d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax) 634d71ae5a4SJacob Faibussowitsch { 6351c3436cfSJed Brown PetscFunctionBegin; 6364782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 637b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmin, 2); 638b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmax, 3); 63909cb0f53SBarry Smith PetscCheck(hmin == (PetscReal)PETSC_CURRENT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin); 64009cb0f53SBarry Smith PetscCheck(hmax == (PetscReal)PETSC_CURRENT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax); 64109cb0f53SBarry Smith if (hmin != (PetscReal)PETSC_CURRENT) adapt->dt_min = hmin; 64209cb0f53SBarry Smith if (hmax != (PetscReal)PETSC_CURRENT) adapt->dt_max = hmax; 643b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 644b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6453c633725SBarry Smith PetscCheck(hmax > hmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum time step %g must greater than minimum time step %g", (double)hmax, (double)hmin); 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6471917a363SLisandro Dalcin } 6481917a363SLisandro Dalcin 6491917a363SLisandro Dalcin /*@ 650bcf0153eSBarry Smith TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the time step controller 6511917a363SLisandro Dalcin 6521917a363SLisandro Dalcin Not Collective 6531917a363SLisandro Dalcin 6544165533cSJose E. Roman Input Parameter: 655bcf0153eSBarry Smith . adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 6561917a363SLisandro Dalcin 6574165533cSJose E. Roman Output Parameters: 6581917a363SLisandro Dalcin + hmin - minimum time step 6591917a363SLisandro Dalcin - hmax - maximum time step 6601917a363SLisandro Dalcin 6611917a363SLisandro Dalcin Level: intermediate 6621917a363SLisandro Dalcin 6631cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 6641917a363SLisandro Dalcin @*/ 665d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax) 666d71ae5a4SJacob Faibussowitsch { 6671917a363SLisandro Dalcin PetscFunctionBegin; 6681917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 6694f572ea9SToby Isaac if (hmin) PetscAssertPointer(hmin, 2); 6704f572ea9SToby Isaac if (hmax) PetscAssertPointer(hmax, 3); 6711917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6721917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6741c3436cfSJed Brown } 6751c3436cfSJed Brown 67614d0ab18SJacob Faibussowitsch /*@C 677bcf0153eSBarry Smith TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options. 67884df9cb4SJed Brown 679c3339decSBarry Smith Collective 68084df9cb4SJed Brown 68114d0ab18SJacob Faibussowitsch Input Parameters: 6822fe279fdSBarry Smith + adapt - the `TSAdapt` context 6832fe279fdSBarry Smith - PetscOptionsObject - object created by `PetscOptionsBegin()` 68484df9cb4SJed Brown 68584df9cb4SJed Brown Options Database Keys: 6861917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 687bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6881917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6891917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6901917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 69123de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 69223de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 69323de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 694de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 695de50f1caSBarry Smith - -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver 69684df9cb4SJed Brown 69784df9cb4SJed Brown Level: advanced 69884df9cb4SJed Brown 699bcf0153eSBarry Smith Note: 700bcf0153eSBarry Smith This function is automatically called by `TSSetFromOptions()` 70184df9cb4SJed Brown 7021cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 703db781477SPatrick Sanan `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 70414d0ab18SJacob Faibussowitsch @*/ 705d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject) 706d71ae5a4SJacob Faibussowitsch { 70784df9cb4SJed Brown char type[256] = TSADAPTBASIC; 70862c23b28SMark PetscReal safety, reject_safety, clip[2], scale, hmin, hmax; 7091c3436cfSJed Brown PetscBool set, flg; 7101917a363SLisandro Dalcin PetscInt two; 71184df9cb4SJed Brown 71284df9cb4SJed Brown PetscFunctionBegin; 713dbbe0bcdSBarry Smith PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 71484df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 7151566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 716d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options"); 7179566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ts_adapt_type", "Algorithm to use for adaptivity", "TSAdaptSetType", TSAdaptList, ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type, type, sizeof(type), &flg)); 71848a46eb9SPierre Jolivet if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type)); 7191917a363SLisandro Dalcin 7209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set)); 7219566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg)); 7221917a363SLisandro Dalcin 7239371c9d4SSatish Balay safety = adapt->safety; 7249371c9d4SSatish Balay reject_safety = adapt->reject_safety; 7259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set)); 7269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg)); 7279566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety)); 7281917a363SLisandro Dalcin 7299371c9d4SSatish Balay two = 2; 7309371c9d4SSatish Balay clip[0] = adapt->clip[0]; 7319371c9d4SSatish Balay clip[1] = adapt->clip[1]; 7329566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set)); 7333c633725SBarry Smith PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip"); 7349566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1])); 7351917a363SLisandro Dalcin 7369371c9d4SSatish Balay hmin = adapt->dt_min; 7379371c9d4SSatish Balay hmax = adapt->dt_max; 7389566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set)); 7399566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg)); 7409566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax)); 7411917a363SLisandro Dalcin 7429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set)); 7439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_glee_use_local", "GLEE adaptor uses local error estimation for step control", "", adapt->glee_use_local, &adapt->glee_use_local, &set)); 744d580f011SEmil Constantinescu 7459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set)); 7469566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale)); 7471917a363SLisandro Dalcin 7489566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL)); 7493c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported"); 7501917a363SLisandro Dalcin 7519566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_adapt_time_step_increase_delay", "Number of timesteps to delay increasing the time step after it has been decreased due to failed solver", "TSAdaptSetTimeStepIncreaseDelay", adapt->timestepjustdecreased_delay, &adapt->timestepjustdecreased_delay, NULL)); 752de50f1caSBarry Smith 7539566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set)); 7549566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt, flg)); 7551917a363SLisandro Dalcin 756dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 757d0609cedSBarry Smith PetscOptionsHeadEnd(); 7583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75984df9cb4SJed Brown } 76084df9cb4SJed Brown 76184df9cb4SJed Brown /*@ 76284df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 76384df9cb4SJed Brown 76420f4b53cSBarry Smith Logically Collective 76584df9cb4SJed Brown 7664165533cSJose E. Roman Input Parameter: 76784df9cb4SJed Brown . adapt - adaptive controller 76884df9cb4SJed Brown 76984df9cb4SJed Brown Level: developer 77084df9cb4SJed Brown 7711cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 77284df9cb4SJed Brown @*/ 773d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 774d71ae5a4SJacob Faibussowitsch { 77584df9cb4SJed Brown PetscFunctionBegin; 7764782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 7779566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates))); 7783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77984df9cb4SJed Brown } 78084df9cb4SJed Brown 78184df9cb4SJed Brown /*@C 78284df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 78384df9cb4SJed Brown 78420f4b53cSBarry Smith Logically Collective; No Fortran Support 78584df9cb4SJed Brown 7864165533cSJose E. Roman Input Parameters: 787bcf0153eSBarry Smith + adapt - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()` 78884df9cb4SJed Brown . name - name of the candidate scheme to add 78984df9cb4SJed Brown . order - order of the candidate scheme 79084df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7918d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 79284df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 79384df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 79484df9cb4SJed Brown 79584df9cb4SJed Brown Level: developer 79684df9cb4SJed Brown 7971cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 79884df9cb4SJed Brown @*/ 799d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse) 800d71ae5a4SJacob Faibussowitsch { 80184df9cb4SJed Brown PetscInt c; 80284df9cb4SJed Brown 80384df9cb4SJed Brown PetscFunctionBegin; 80484df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 80563a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order); 80684df9cb4SJed Brown if (inuse) { 8073c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 80884df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 80984df9cb4SJed Brown } 8101c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 8111c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 812bbd56ea5SKarl Rupp 81384df9cb4SJed Brown adapt->candidates.name[c] = name; 81484df9cb4SJed Brown adapt->candidates.order[c] = order; 81584df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 8168d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 81784df9cb4SJed Brown adapt->candidates.cost[c] = cost; 81884df9cb4SJed Brown adapt->candidates.n++; 8193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82084df9cb4SJed Brown } 82184df9cb4SJed Brown 8228d59e960SJed Brown /*@C 8238d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8248d59e960SJed Brown 8258d59e960SJed Brown Not Collective 8268d59e960SJed Brown 8274165533cSJose E. Roman Input Parameter: 8288d59e960SJed Brown . adapt - time step adaptivity context 8298d59e960SJed Brown 8304165533cSJose E. Roman Output Parameters: 8318d59e960SJed Brown + n - number of candidate schemes, always at least 1 8328d59e960SJed Brown . order - the order of each candidate scheme 8338d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8348d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8358d59e960SJed Brown - cost - the relative cost of each scheme 8368d59e960SJed Brown 8378d59e960SJed Brown Level: developer 8388d59e960SJed Brown 8398d59e960SJed Brown Note: 8408d59e960SJed Brown The current scheme is always returned in the first slot 8418d59e960SJed Brown 8421cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 8438d59e960SJed Brown @*/ 844d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost) 845d71ae5a4SJacob Faibussowitsch { 8468d59e960SJed Brown PetscFunctionBegin; 8478d59e960SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 8488d59e960SJed Brown if (n) *n = adapt->candidates.n; 8498d59e960SJed Brown if (order) *order = adapt->candidates.order; 8508d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8518d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8528d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8548d59e960SJed Brown } 8558d59e960SJed Brown 85684df9cb4SJed Brown /*@C 85784df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 85884df9cb4SJed Brown 859c3339decSBarry Smith Collective 86084df9cb4SJed Brown 8614165533cSJose E. Roman Input Parameters: 862da81f932SPierre Jolivet + adapt - adaptive controller 86397bb3fdcSJose E. Roman . ts - time stepper 86484df9cb4SJed Brown - h - current step size 86584df9cb4SJed Brown 8664165533cSJose E. Roman Output Parameters: 8671566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 86884df9cb4SJed Brown . next_h - step size to use for the next step 869bcf0153eSBarry Smith - accept - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size 8701c3436cfSJed Brown 87184df9cb4SJed Brown Level: developer 87284df9cb4SJed Brown 873bcf0153eSBarry Smith Note: 874bcf0153eSBarry Smith The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is 875bcf0153eSBarry Smith being retried after an initial rejection. 876bcf0153eSBarry Smith 8771cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 87884df9cb4SJed Brown @*/ 879d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept) 880d71ae5a4SJacob Faibussowitsch { 8811566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8821566a47fSLisandro Dalcin PetscInt scheme = 0; 8830b99f514SJed Brown PetscReal wlte = -1.0; 8845e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8855e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 88684df9cb4SJed Brown 88784df9cb4SJed Brown PetscFunctionBegin; 88884df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 88984df9cb4SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 8904f572ea9SToby Isaac if (next_sc) PetscAssertPointer(next_sc, 4); 8914f572ea9SToby Isaac PetscAssertPointer(next_h, 5); 8924f572ea9SToby Isaac PetscAssertPointer(accept, 6); 8931566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8941566a47fSLisandro Dalcin 8951566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events */ 896ca4445c7SIlya Fursov if (ts->event && ts->event->processing) { 8971566a47fSLisandro Dalcin *next_h = h; 8981566a47fSLisandro Dalcin *accept = PETSC_TRUE; 899fe4ad979SIlya Fursov if (adapt->monitor) { 900fe4ad979SIlya Fursov PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 901fe4ad979SIlya Fursov 902fe4ad979SIlya Fursov if (ts->event->iterctr == 0) { 903fe4ad979SIlya Fursov /* 904fe4ad979SIlya Fursov An event has been found, now finalising the event processing: performing the 1st and 2nd post-event steps. 905fe4ad979SIlya Fursov Entering this if-branch means both these steps (set to either PETSC_DECIDE or numerical value) are managed 906fe4ad979SIlya Fursov by the event handler. In this case the 1st post-event step is always accepted, without interference of TSAdapt. 907fe4ad979SIlya Fursov Note: if the 2nd post-event step is not managed by the event handler (e.g. given 1st = numerical, 2nd = PETSC_DECIDE), 908fe4ad979SIlya Fursov this if-branch is not entered, and TSAdapt may reject/adjust the proposed 1st post-event step. 909fe4ad979SIlya Fursov */ 910fe4ad979SIlya Fursov PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt does not interfere, step %3" PetscInt_FMT " accepted. Processing post-event steps: 1-st accepted just now, 2-nd yet to come\n", ts->steps)); 911fe4ad979SIlya Fursov } else PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt does not interfere, step %3" PetscInt_FMT " accepted. Event handling in progress\n", ts->steps)); 912fe4ad979SIlya Fursov 913fe4ad979SIlya Fursov PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 914fe4ad979SIlya Fursov } 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9161566a47fSLisandro Dalcin } 9171566a47fSLisandro Dalcin 918dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter); 91963a3b9bcSJacob Faibussowitsch PetscCheck(scheme >= 0 && (ncandidates <= 0 || scheme < ncandidates), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Chosen scheme %" PetscInt_FMT " not in valid range 0..%" PetscInt_FMT, scheme, ncandidates - 1); 9203c633725SBarry Smith PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h); 9211566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 9221566a47fSLisandro Dalcin 92349354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 92436b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 925ca4445c7SIlya Fursov PetscReal t = ts->ptime + ts->time_step, tend, tmax, h1, hmax; 926fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 927fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 928ca4445c7SIlya Fursov PetscReal eps = 10 * PETSC_MACHINE_EPSILON; 9294a658b32SHong Zhang 930fe4ad979SIlya Fursov /* 931fe4ad979SIlya Fursov Logic in using 'dt_span_cached': 932fe4ad979SIlya Fursov 1. It always overrides *next_h, except (any of): 933fe4ad979SIlya Fursov a) the current step was rejected, 934fe4ad979SIlya Fursov b) the adaptor proposed to decrease the next step, 935fe4ad979SIlya Fursov c) the adaptor proposed *next_h > dt_span_cached. 936*136cf249SJames Wright 2. If *next_h was adjusted by eval_times points (or the final point): 937fe4ad979SIlya Fursov -- when dt_span_cached is filled (>0), it keeps its value, 938fe4ad979SIlya Fursov -- when dt_span_cached is clear (==0), it gets the unadjusted version of *next_h. 939fe4ad979SIlya Fursov 3. If *next_h was not adjusted as in (2), dt_span_cached is cleared. 940fe4ad979SIlya Fursov Note, if a combination (1.b || 1.c) && (3) takes place, this means that 941fe4ad979SIlya Fursov dt_span_cached remains unused at the moment of clearing. 942fe4ad979SIlya Fursov If (1.a) takes place, dt_span_cached keeps its value. 943fe4ad979SIlya Fursov Also, dt_span_cached can be updated by the event handler, see tsevent.c. 944fe4ad979SIlya Fursov */ 945*136cf249SJames Wright if (h <= *next_h && *next_h <= adapt->dt_eval_times_cached) *next_h = adapt->dt_eval_times_cached; /* try employing the cache */ 946ca4445c7SIlya Fursov h1 = *next_h; 947ca4445c7SIlya Fursov tend = t + h1; 948ca4445c7SIlya Fursov 949*136cf249SJames Wright if (ts->eval_times && ts->eval_times->time_point_idx < ts->eval_times->num_time_points) { 950*136cf249SJames Wright PetscCheck(ts->eval_times->worktol == 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_PLIB, "Unexpected state (tspan->worktol != 0) in TSAdaptChoose()"); 951*136cf249SJames Wright ts->eval_times->worktol = ts->eval_times->reltol * h1 + ts->eval_times->abstol; 952*136cf249SJames Wright if (PetscIsCloseAtTol(t, ts->eval_times->time_points[ts->eval_times->time_point_idx], ts->eval_times->worktol, 0)) /* hit a span time point */ 953*136cf249SJames Wright if (ts->eval_times->time_point_idx + 1 < ts->eval_times->num_time_points) tmax = ts->eval_times->time_points[ts->eval_times->time_point_idx + 1]; 9544a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */ 955*136cf249SJames Wright else tmax = ts->eval_times->time_points[ts->eval_times->time_point_idx]; 9564a658b32SHong Zhang } else tmax = ts->max_time; 957ca4445c7SIlya Fursov tmax = PetscMin(tmax, ts->max_time); 9584a658b32SHong Zhang hmax = tmax - t; 959ca4445c7SIlya Fursov PetscCheck((hmax > eps) || (PetscAbsReal(hmax) <= eps && PetscIsCloseAtTol(t, ts->max_time, eps, 0)), PetscObjectComm((PetscObject)adapt), PETSC_ERR_PLIB, "Unexpected state: bad hmax in TSAdaptChoose()"); 960ca4445c7SIlya Fursov 96136b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 962ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * b > hmax) *next_h = hmax / 2; 963ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * a > hmax) *next_h = hmax; 964*136cf249SJames Wright if (ts->eval_times && h1 != *next_h && !adapt->dt_eval_times_cached) adapt->dt_eval_times_cached = h1; /* cache the step size if it is to be changed */ 965*136cf249SJames Wright if (ts->eval_times && h1 == *next_h && adapt->dt_eval_times_cached) adapt->dt_eval_times_cached = 0; /* clear the cache if the step size is unchanged */ 966e1db57b0SHong Zhang } 9671c3436cfSJed Brown if (adapt->monitor) { 9681566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 9699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9700b99f514SJed Brown if (wlte < 0) { 9719371c9d4SSatish Balay PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s %s %" PetscInt_FMT ":%s step %3" PetscInt_FMT " %s t=%-11g+%10.3e dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ((PetscObject)ts)->type_name, scheme, sc_name, ts->steps, *accept ? "accepted" : "rejected", 9729371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h)); 9730b99f514SJed Brown } else { 9749371c9d4SSatish Balay PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s %s %" PetscInt_FMT ":%s step %3" PetscInt_FMT " %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g wltea=%5.3g wlter=%5.3g\n", ((PetscObject)adapt)->type_name, ((PetscObject)ts)->type_name, scheme, sc_name, ts->steps, *accept ? "accepted" : "rejected", 9759371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter)); 9760b99f514SJed Brown } 9779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9781c3436cfSJed Brown } 9793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98084df9cb4SJed Brown } 98184df9cb4SJed Brown 98297335746SJed Brown /*@ 983de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 984de50f1caSBarry Smith before increasing the time step. 985de50f1caSBarry Smith 986c3339decSBarry Smith Logicially Collective 987de50f1caSBarry Smith 9884165533cSJose E. Roman Input Parameters: 989de50f1caSBarry Smith + adapt - adaptive controller context 990de50f1caSBarry Smith - cnt - the number of timesteps 991de50f1caSBarry Smith 992de50f1caSBarry Smith Options Database Key: 993de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 994de50f1caSBarry Smith 995de50f1caSBarry Smith Level: advanced 996de50f1caSBarry Smith 997bcf0153eSBarry Smith Notes: 998bcf0153eSBarry Smith This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 999bcf0153eSBarry Smith 1000bcf0153eSBarry Smith The successful use of this option is problem dependent 1001bcf0153eSBarry Smith 1002b43aa488SJacob Faibussowitsch Developer Notes: 1003bcf0153eSBarry Smith There is no theory to support this option 1004bcf0153eSBarry Smith 10051cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt` 1006de50f1caSBarry Smith @*/ 1007d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt) 1008d71ae5a4SJacob Faibussowitsch { 1009de50f1caSBarry Smith PetscFunctionBegin; 1010de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 10113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1012de50f1caSBarry Smith } 1013de50f1caSBarry Smith 1014de50f1caSBarry Smith /*@ 10156bc98fa9SBarry Smith TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails or solution vector is infeasible) 101697335746SJed Brown 1017c3339decSBarry Smith Collective 101897335746SJed Brown 10194165533cSJose E. Roman Input Parameters: 102097335746SJed Brown + adapt - adaptive controller context 1021b295832fSPierre Barbier de Reuille . ts - time stepper 1022b295832fSPierre Barbier de Reuille . t - Current simulation time 1023b295832fSPierre Barbier de Reuille - Y - Current solution vector 102497335746SJed Brown 10254165533cSJose E. Roman Output Parameter: 1026bcf0153eSBarry Smith . accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject 102797335746SJed Brown 102897335746SJed Brown Level: developer 102997335746SJed Brown 10301cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt` 103197335746SJed Brown @*/ 1032d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) 1033d71ae5a4SJacob Faibussowitsch { 10341566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 10356c6709e3SStefano Zampini PetscBool func_accept, snes_div_func; 103697335746SJed Brown 103797335746SJed Brown PetscFunctionBegin; 10384782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 10394782b174SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 10404f572ea9SToby Isaac PetscAssertPointer(accept, 5); 10411566a47fSLisandro Dalcin 10426c6709e3SStefano Zampini PetscCall(TSFunctionDomainError(ts, t, Y, &func_accept)); 10439566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 10446c6709e3SStefano Zampini snes_div_func = (PetscBool)(snesreason == SNES_DIVERGED_FUNCTION_DOMAIN); 10456c6709e3SStefano Zampini if (func_accept && snesreason < 0 && !snes_div_func) { 104697335746SJed Brown *accept = PETSC_FALSE; 10476c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failure: %s\n", ts->steps, SNESConvergedReasons[snesreason])); 104809cb0f53SBarry Smith if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures != PETSC_UNLIMITED) { 104997335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 105063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 105197335746SJed Brown if (adapt->monitor) { 10529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10539371c9d4SSatish Balay PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected t=%-11g+%10.3e, nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed\n", ((PetscObject)adapt)->type_name, ts->steps, 10549371c9d4SSatish Balay (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures)); 10559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 105697335746SJed Brown } 1057cb9d8021SPierre Barbier de Reuille } 1058cb9d8021SPierre Barbier de Reuille } else { 10596c6709e3SStefano Zampini *accept = (PetscBool)(func_accept && !snes_div_func); 10606c6709e3SStefano Zampini if (*accept && adapt->checkstage) PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept)); 10616bc98fa9SBarry Smith if (!*accept) { 10626c6709e3SStefano Zampini const char *user_func = !func_accept ? "TSSetFunctionDomainError()" : "TSAdaptSetCheckStage"; 10636c6709e3SStefano Zampini const char *snes_err = "SNES invalid function domain"; 10646c6709e3SStefano Zampini const char *err_msg = snes_div_func && func_accept ? snes_err : user_func; 10656c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by %s\n", ts->steps, err_msg)); 10666bc98fa9SBarry Smith if (adapt->monitor) { 10679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10686c6709e3SStefano Zampini PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected by %s\n", ((PetscObject)adapt)->type_name, ts->steps, err_msg)); 10699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10706bc98fa9SBarry Smith } 10716bc98fa9SBarry Smith } 1072cb9d8021SPierre Barbier de Reuille } 1073cb9d8021SPierre Barbier de Reuille 107457508eceSPierre Jolivet if (!*accept && !ts->reason) { 10751566a47fSLisandro Dalcin PetscReal dt, new_dt; 10769566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1077cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 10789566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, new_dt)); 1079de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 108097335746SJed Brown if (adapt->monitor) { 10819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10826c6709e3SStefano Zampini PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected (SNES reason %s) t=%-11g+%10.3e retrying with dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ts->steps, SNESConvergedReasons[snesreason], 10836c6709e3SStefano Zampini (double)ts->ptime, (double)dt, (double)new_dt)); 10849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 108597335746SJed Brown } 108697335746SJed Brown } 10873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 108897335746SJed Brown } 108997335746SJed Brown 109084df9cb4SJed Brown /*@ 109184df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 109284df9cb4SJed Brown 1093d083f849SBarry Smith Collective 109484df9cb4SJed Brown 109584df9cb4SJed Brown Input Parameter: 109684df9cb4SJed Brown . comm - The communicator 109784df9cb4SJed Brown 109884df9cb4SJed Brown Output Parameter: 1099b43aa488SJacob Faibussowitsch . inadapt - new `TSAdapt` object 110084df9cb4SJed Brown 110184df9cb4SJed Brown Level: developer 110284df9cb4SJed Brown 1103bcf0153eSBarry Smith Note: 1104bcf0153eSBarry Smith `TSAdapt` creation is handled by `TS`, so users should not need to call this function. 110584df9cb4SJed Brown 11061cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 110784df9cb4SJed Brown @*/ 1108d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt) 1109d71ae5a4SJacob Faibussowitsch { 111084df9cb4SJed Brown TSAdapt adapt; 111184df9cb4SJed Brown 111284df9cb4SJed Brown PetscFunctionBegin; 11134f572ea9SToby Isaac PetscAssertPointer(inadapt, 2); 11149566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 11153b3bcf4cSLisandro Dalcin 11169566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView)); 1117bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 11181917a363SLisandro Dalcin adapt->safety = 0.9; 11191917a363SLisandro Dalcin adapt->reject_safety = 0.5; 11201917a363SLisandro Dalcin adapt->clip[0] = 0.1; 11211917a363SLisandro Dalcin adapt->clip[1] = 10.; 11221c3436cfSJed Brown adapt->dt_min = 1e-20; 11231917a363SLisandro Dalcin adapt->dt_max = 1e+20; 11241c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1125d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 112697335746SJed Brown adapt->scale_solve_failed = 0.25; 1127fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1128fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1129fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1130fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 11317619abb3SShri adapt->wnormtype = NORM_2; 1132de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 113384df9cb4SJed Brown *inadapt = adapt; 11343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113584df9cb4SJed Brown } 1136