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 36*b44f4de4SBarry Smith .vb 37*b44f4de4SBarry Smith TSAdaptSetType(ts, "my_scheme") 38*b44f4de4SBarry Smith .ve 391c84c290SBarry Smith or at runtime via the option 40*b44f4de4SBarry Smith .vb 41*b44f4de4SBarry Smith -ts_adapt_type my_scheme 42*b44f4de4SBarry Smith .ve 4384df9cb4SJed Brown 441cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterAll()` 4584df9cb4SJed Brown @*/ 46d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt)) 47d71ae5a4SJacob Faibussowitsch { 4884df9cb4SJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 509566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function)); 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5284df9cb4SJed Brown } 5384df9cb4SJed Brown 5484df9cb4SJed Brown /*@C 552fe279fdSBarry Smith TSAdaptRegisterAll - Registers all of the adaptivity schemes in `TSAdapt` 5684df9cb4SJed Brown 5784df9cb4SJed Brown Not Collective 5884df9cb4SJed Brown 5984df9cb4SJed Brown Level: advanced 6084df9cb4SJed Brown 611cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterDestroy()` 6284df9cb4SJed Brown @*/ 63d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegisterAll(void) 64d71ae5a4SJacob Faibussowitsch { 6584df9cb4SJed Brown PetscFunctionBegin; 663ba16761SJacob Faibussowitsch if (TSAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 670f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 689566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None)); 699566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic)); 709566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP)); 719566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL)); 729566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE)); 739566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History)); 743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7584df9cb4SJed Brown } 7684df9cb4SJed Brown 7784df9cb4SJed Brown /*@C 782fe279fdSBarry Smith TSAdaptFinalizePackage - This function destroys everything in the `TS` package. It is 792fe279fdSBarry Smith called from `PetscFinalize()`. 8084df9cb4SJed Brown 8184df9cb4SJed Brown Level: developer 8284df9cb4SJed Brown 831cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()` 8484df9cb4SJed Brown @*/ 85d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptFinalizePackage(void) 86d71ae5a4SJacob Faibussowitsch { 8784df9cb4SJed Brown PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSAdaptList)); 8984df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9084df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9284df9cb4SJed Brown } 9384df9cb4SJed Brown 9484df9cb4SJed Brown /*@C 952fe279fdSBarry Smith TSAdaptInitializePackage - This function initializes everything in the `TSAdapt` package. It is 962fe279fdSBarry Smith called from `TSInitializePackage()`. 9784df9cb4SJed Brown 9884df9cb4SJed Brown Level: developer 9984df9cb4SJed Brown 1001cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()` 10184df9cb4SJed Brown @*/ 102d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptInitializePackage(void) 103d71ae5a4SJacob Faibussowitsch { 10484df9cb4SJed Brown PetscFunctionBegin; 1053ba16761SJacob Faibussowitsch if (TSAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 10684df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 1079566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID)); 1089566063dSJacob Faibussowitsch PetscCall(TSAdaptRegisterAll()); 1099566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage)); 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11184df9cb4SJed Brown } 11284df9cb4SJed Brown 113cc4c1da9SBarry Smith /*@ 114bcf0153eSBarry Smith TSAdaptSetType - sets the approach used for the error adapter 1157eef6055SBarry Smith 11620f4b53cSBarry Smith Logicially Collective 1177eef6055SBarry Smith 118d8d19677SJose E. Roman Input Parameters: 11958b119d1SBarry Smith + adapt - the `TS` adapter, most likely obtained with `TSGetAdapt()` 120bcf0153eSBarry Smith - type - one of the `TSAdaptType` 1217eef6055SBarry Smith 122bcf0153eSBarry Smith Options Database Key: 123ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1247eef6055SBarry Smith 1257eef6055SBarry Smith Level: intermediate 1267eef6055SBarry Smith 127b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()` 1287eef6055SBarry Smith @*/ 129d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type) 130d71ae5a4SJacob Faibussowitsch { 131ef749922SLisandro Dalcin PetscBool match; 1325f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSAdapt); 13384df9cb4SJed Brown 13484df9cb4SJed Brown PetscFunctionBegin; 1354782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1364f572ea9SToby Isaac PetscAssertPointer(type, 2); 1379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match)); 1383ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1399566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSAdaptList, type, &r)); 1406adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type); 141dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, destroy); 1429566063dSJacob Faibussowitsch PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps))); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 1449566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14684df9cb4SJed Brown } 14784df9cb4SJed Brown 148cc4c1da9SBarry Smith /*@ 149a1cb98faSBarry Smith TSAdaptGetType - gets the `TS` adapter method type (as a string). 150d0288e4fSLisandro Dalcin 151d0288e4fSLisandro Dalcin Not Collective 152d0288e4fSLisandro Dalcin 153d0288e4fSLisandro Dalcin Input Parameter: 154a1cb98faSBarry Smith . adapt - The `TS` adapter, most likely obtained with `TSGetAdapt()` 155d0288e4fSLisandro Dalcin 156d0288e4fSLisandro Dalcin Output Parameter: 157a1cb98faSBarry Smith . type - The name of `TS` adapter method 158d0288e4fSLisandro Dalcin 159d0288e4fSLisandro Dalcin Level: intermediate 160d0288e4fSLisandro Dalcin 161a1cb98faSBarry Smith .seealso: `TSAdapt`, `TSAdaptType`, `TSAdaptSetType()` 162d0288e4fSLisandro Dalcin @*/ 163d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type) 164d71ae5a4SJacob Faibussowitsch { 165d0288e4fSLisandro Dalcin PetscFunctionBegin; 166d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1674f572ea9SToby Isaac PetscAssertPointer(type, 2); 168d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 1693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 170d0288e4fSLisandro Dalcin } 171d0288e4fSLisandro Dalcin 172d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[]) 173d71ae5a4SJacob Faibussowitsch { 17484df9cb4SJed Brown PetscFunctionBegin; 1754782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17884df9cb4SJed Brown } 17984df9cb4SJed Brown 180ffeef943SBarry Smith /*@ 1812fe279fdSBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with `TSAdaptView()`. 182ad6bc421SBarry Smith 183c3339decSBarry Smith Collective 184ad6bc421SBarry Smith 185ad6bc421SBarry Smith Input Parameters: 186b43aa488SJacob Faibussowitsch + adapt - the newly loaded `TSAdapt`, this needs to have been created with `TSAdaptCreate()` or 187bcf0153eSBarry Smith some related function before a call to `TSAdaptLoad()`. 188bcf0153eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or 189bcf0153eSBarry Smith HDF5 file viewer, obtained from `PetscViewerHDF5Open()` 190ad6bc421SBarry Smith 191ad6bc421SBarry Smith Level: intermediate 192ad6bc421SBarry Smith 193bcf0153eSBarry Smith Note: 194bcf0153eSBarry Smith The type is determined by the data in the file, any type set into the `TSAdapt` before this call is ignored. 195ad6bc421SBarry Smith 1961cc06b55SBarry Smith .seealso: [](ch_ts), `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()`, `TSAdapt` 197ad6bc421SBarry Smith @*/ 198d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer) 199d71ae5a4SJacob Faibussowitsch { 200ad6bc421SBarry Smith PetscBool isbinary; 201ad6bc421SBarry Smith char type[256]; 202ad6bc421SBarry Smith 203ad6bc421SBarry Smith PetscFunctionBegin; 2044782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 205ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 2073c633725SBarry Smith PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 208ad6bc421SBarry Smith 2099566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 2109566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, type)); 211dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, load, viewer); 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 213ad6bc421SBarry Smith } 214ad6bc421SBarry Smith 215d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer) 216d71ae5a4SJacob Faibussowitsch { 2171c167fc2SEmil Constantinescu PetscBool iascii, isbinary, isnone, isglee; 21884df9cb4SJed Brown 21984df9cb4SJed Brown PetscFunctionBegin; 2204782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 2219566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer)); 2224782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2234782b174SLisandro Dalcin PetscCheckSameComm(adapt, 1, viewer, 2); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 22684df9cb4SJed Brown if (iascii) { 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone)); 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee)); 2301917a363SLisandro Dalcin if (!isnone) { 2319566063dSJacob Faibussowitsch if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, " always accepting steps\n")); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " safety factor %g\n", (double)adapt->safety)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " extra safety factor after step rejection %g\n", (double)adapt->reject_safety)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest increase %g\n", (double)adapt->clip[1])); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest decrease %g\n", (double)adapt->clip[0])); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum allowed timestep %g\n", (double)adapt->dt_max)); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " minimum allowed timestep %g\n", (double)adapt->dt_min)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max)); 2391c167fc2SEmil Constantinescu } 2401c167fc2SEmil Constantinescu if (isglee) { 2411c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses local error control\n")); 2431c167fc2SEmil Constantinescu } else { 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses global error control\n")); 2451c167fc2SEmil Constantinescu } 2461917a363SLisandro Dalcin } 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 248dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, view, viewer); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 250ad6bc421SBarry Smith } else if (isbinary) { 251ad6bc421SBarry Smith char type[256]; 252ad6bc421SBarry Smith 253ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 2549566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256)); 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 256dbbe0bcdSBarry Smith } else PetscTryTypeMethod(adapt, view, viewer); 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25884df9cb4SJed Brown } 25984df9cb4SJed Brown 260099af0a3SLisandro Dalcin /*@ 261bcf0153eSBarry Smith TSAdaptReset - Resets a `TSAdapt` context to its defaults 262099af0a3SLisandro Dalcin 263c3339decSBarry Smith Collective 264099af0a3SLisandro Dalcin 265099af0a3SLisandro Dalcin Input Parameter: 266bcf0153eSBarry Smith . adapt - the `TSAdapt` context obtained from `TSGetAdapt()` or `TSAdaptCreate()` 267099af0a3SLisandro Dalcin 268099af0a3SLisandro Dalcin Level: developer 269099af0a3SLisandro Dalcin 2701cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdapt`, `TSAdaptCreate()`, `TSAdaptDestroy()` 271099af0a3SLisandro Dalcin @*/ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptReset(TSAdapt adapt) 273d71ae5a4SJacob Faibussowitsch { 274099af0a3SLisandro Dalcin PetscFunctionBegin; 275099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 276dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, reset); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278099af0a3SLisandro Dalcin } 279099af0a3SLisandro Dalcin 280d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 281d71ae5a4SJacob Faibussowitsch { 28284df9cb4SJed Brown PetscFunctionBegin; 2833ba16761SJacob Faibussowitsch if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS); 28484df9cb4SJed Brown PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1); 285f4f49eeaSPierre Jolivet if (--((PetscObject)*adapt)->refct > 0) { 2869371c9d4SSatish Balay *adapt = NULL; 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2889371c9d4SSatish Balay } 289099af0a3SLisandro Dalcin 2909566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt)); 291099af0a3SLisandro Dalcin 292f4f49eeaSPierre Jolivet PetscTryTypeMethod(*adapt, destroy); 2939566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 2949566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29684df9cb4SJed Brown } 29784df9cb4SJed Brown 2981c3436cfSJed Brown /*@ 2991c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3001c3436cfSJed Brown 301c3339decSBarry Smith Collective 3021c3436cfSJed Brown 3034165533cSJose E. Roman Input Parameters: 3041c3436cfSJed Brown + adapt - adaptive controller context 305bcf0153eSBarry Smith - flg - `PETSC_TRUE` to active a monitor, `PETSC_FALSE` to disable 3061c3436cfSJed Brown 307bcf0153eSBarry Smith Options Database Key: 308ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 309bf997491SLisandro Dalcin 3101c3436cfSJed Brown Level: intermediate 3111c3436cfSJed Brown 3121cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 3131c3436cfSJed Brown @*/ 314d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg) 315d71ae5a4SJacob Faibussowitsch { 3161c3436cfSJed Brown PetscFunctionBegin; 3174782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3184782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flg, 2); 3191c3436cfSJed Brown if (flg) { 3209566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor)); 3211c3436cfSJed Brown } else { 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor)); 3231c3436cfSJed Brown } 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3251c3436cfSJed Brown } 3261c3436cfSJed Brown 3270873808bSJed Brown /*@C 328bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3290873808bSJed Brown 3302fe279fdSBarry Smith Logically Collective 3310873808bSJed Brown 3324165533cSJose E. Roman Input Parameters: 3330873808bSJed Brown + adapt - adaptive controller context 3340873808bSJed Brown - func - stage check function 3350873808bSJed Brown 336b43aa488SJacob Faibussowitsch Calling sequence: 3370873808bSJed Brown + adapt - adaptive controller context 3380873808bSJed Brown . ts - time stepping context 33914d0ab18SJacob Faibussowitsch . t - current time 34014d0ab18SJacob Faibussowitsch . Y - current solution vector 3410873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3420873808bSJed Brown 3430873808bSJed Brown Level: advanced 3440873808bSJed Brown 3451cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 3460873808bSJed Brown @*/ 34714d0ab18SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)) 348d71ae5a4SJacob Faibussowitsch { 3490873808bSJed Brown PetscFunctionBegin; 3500873808bSJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 35168ae941aSLisandro Dalcin adapt->checkstage = func; 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3530873808bSJed Brown } 3540873808bSJed Brown 3551c3436cfSJed Brown /*@ 356bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 357bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 358bf997491SLisandro Dalcin 3592fe279fdSBarry Smith Logically Collective 360bf997491SLisandro Dalcin 3614165533cSJose E. Roman Input Parameters: 362bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 363bf997491SLisandro Dalcin - flag - whether to always accept steps 364bf997491SLisandro Dalcin 365bcf0153eSBarry Smith Options Database Key: 366ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 367bf997491SLisandro Dalcin 368bf997491SLisandro Dalcin Level: intermediate 369bf997491SLisandro Dalcin 3701cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()` 371bf997491SLisandro Dalcin @*/ 372d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag) 373d71ae5a4SJacob Faibussowitsch { 374bf997491SLisandro Dalcin PetscFunctionBegin; 375bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 376bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flag, 2); 377bf997491SLisandro Dalcin adapt->always_accept = flag; 3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 379bf997491SLisandro Dalcin } 380bf997491SLisandro Dalcin 381bf997491SLisandro Dalcin /*@ 382bcf0153eSBarry Smith TSAdaptSetSafety - Set safety factors for time step adaptor 3831c3436cfSJed Brown 3842fe279fdSBarry Smith Logically Collective 3851917a363SLisandro Dalcin 3864165533cSJose E. Roman Input Parameters: 3871917a363SLisandro Dalcin + adapt - adaptive controller context 3881917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 389ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 3901917a363SLisandro Dalcin 3911917a363SLisandro Dalcin Options Database Keys: 392ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 393ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 3941917a363SLisandro Dalcin 3951917a363SLisandro Dalcin Level: intermediate 3961917a363SLisandro Dalcin 39709cb0f53SBarry Smith Note: 39809cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 39909cb0f53SBarry Smith 40009cb0f53SBarry Smith Fortran Note: 40109cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 40209cb0f53SBarry Smith 4031cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 4041917a363SLisandro Dalcin @*/ 405d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety) 406d71ae5a4SJacob Faibussowitsch { 4071917a363SLisandro Dalcin PetscFunctionBegin; 4081917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4091917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, safety, 2); 4101917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, reject_safety, 3); 41109cb0f53SBarry Smith PetscCheck(safety == (PetscReal)PETSC_CURRENT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety); 41209cb0f53SBarry 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); 41309cb0f53SBarry 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); 41409cb0f53SBarry 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); 41509cb0f53SBarry Smith if (safety != (PetscReal)PETSC_CURRENT) adapt->safety = safety; 41609cb0f53SBarry Smith if (reject_safety != (PetscReal)PETSC_CURRENT) adapt->reject_safety = reject_safety; 4173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4181917a363SLisandro Dalcin } 4191917a363SLisandro Dalcin 4201917a363SLisandro Dalcin /*@ 421bcf0153eSBarry Smith TSAdaptGetSafety - Get safety factors for time step adapter 4221917a363SLisandro Dalcin 4231917a363SLisandro Dalcin Not Collective 4241917a363SLisandro Dalcin 4254165533cSJose E. Roman Input Parameter: 4261917a363SLisandro Dalcin . adapt - adaptive controller context 4271917a363SLisandro Dalcin 4284165533cSJose E. Roman Output Parameters: 429b43aa488SJacob Faibussowitsch + safety - safety factor relative to target error/stability goal 430b43aa488SJacob Faibussowitsch - reject_safety - extra safety factor to apply if the last step was rejected 4311917a363SLisandro Dalcin 4321917a363SLisandro Dalcin Level: intermediate 4331917a363SLisandro Dalcin 4341cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 4351917a363SLisandro Dalcin @*/ 436d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety) 437d71ae5a4SJacob Faibussowitsch { 4381917a363SLisandro Dalcin PetscFunctionBegin; 4391917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4404f572ea9SToby Isaac if (safety) PetscAssertPointer(safety, 2); 4414f572ea9SToby Isaac if (reject_safety) PetscAssertPointer(reject_safety, 3); 4421917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4431917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4451917a363SLisandro Dalcin } 4461917a363SLisandro Dalcin 4471917a363SLisandro Dalcin /*@ 448bcf0153eSBarry Smith TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms 449bcf0153eSBarry Smith for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components. 45076cddca1SEmil Constantinescu 4512fe279fdSBarry Smith Logically Collective 45276cddca1SEmil Constantinescu 4534165533cSJose E. Roman Input Parameters: 45476cddca1SEmil Constantinescu + adapt - adaptive controller context 45576cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 45676cddca1SEmil Constantinescu 457bcf0153eSBarry Smith Options Database Key: 45876cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 45976cddca1SEmil Constantinescu 46076cddca1SEmil Constantinescu Level: intermediate 46176cddca1SEmil Constantinescu 4621cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 46376cddca1SEmil Constantinescu @*/ 464d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore) 465d71ae5a4SJacob Faibussowitsch { 46676cddca1SEmil Constantinescu PetscFunctionBegin; 46776cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 46876cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt, max_ignore, 2); 46976cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47176cddca1SEmil Constantinescu } 47276cddca1SEmil Constantinescu 47376cddca1SEmil Constantinescu /*@ 474bcf0153eSBarry Smith TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms 475bcf0153eSBarry Smith for time step adaptivity (in absolute value). 47676cddca1SEmil Constantinescu 47776cddca1SEmil Constantinescu Not Collective 47876cddca1SEmil Constantinescu 4794165533cSJose E. Roman Input Parameter: 48076cddca1SEmil Constantinescu . adapt - adaptive controller context 48176cddca1SEmil Constantinescu 4824165533cSJose E. Roman Output Parameter: 48376cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 48476cddca1SEmil Constantinescu 48576cddca1SEmil Constantinescu Level: intermediate 48676cddca1SEmil Constantinescu 4871cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 48876cddca1SEmil Constantinescu @*/ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore) 490d71ae5a4SJacob Faibussowitsch { 49176cddca1SEmil Constantinescu PetscFunctionBegin; 49276cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4934f572ea9SToby Isaac PetscAssertPointer(max_ignore, 2); 49476cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 4953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49676cddca1SEmil Constantinescu } 49776cddca1SEmil Constantinescu 49876cddca1SEmil Constantinescu /*@ 499bcf0153eSBarry Smith TSAdaptSetClip - Sets the admissible decrease/increase factor in step size in the time step adapter 5001917a363SLisandro Dalcin 501c3339decSBarry Smith Logically collective 5021917a363SLisandro Dalcin 5034165533cSJose E. Roman Input Parameters: 5041917a363SLisandro Dalcin + adapt - adaptive controller context 5051917a363SLisandro Dalcin . low - admissible decrease factor 506ec18b7bcSLisandro Dalcin - high - admissible increase factor 5071917a363SLisandro Dalcin 508bcf0153eSBarry Smith Options Database Key: 509ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5101917a363SLisandro Dalcin 5111917a363SLisandro Dalcin Level: intermediate 5121917a363SLisandro Dalcin 51309cb0f53SBarry Smith Note: 51409cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 51509cb0f53SBarry Smith 51609cb0f53SBarry Smith Fortran Note: 51709cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 51809cb0f53SBarry Smith 5191cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 5201917a363SLisandro Dalcin @*/ 521d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high) 522d71ae5a4SJacob Faibussowitsch { 5231917a363SLisandro Dalcin PetscFunctionBegin; 5241917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 5251917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, low, 2); 5261917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, high, 3); 52709cb0f53SBarry Smith PetscCheck(low == (PetscReal)PETSC_CURRENT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low); 52809cb0f53SBarry 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); 52909cb0f53SBarry 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); 53009cb0f53SBarry Smith if (low != (PetscReal)PETSC_CURRENT) adapt->clip[0] = low; 53109cb0f53SBarry Smith if (high != (PetscReal)PETSC_CURRENT) adapt->clip[1] = high; 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5331917a363SLisandro Dalcin } 5341917a363SLisandro Dalcin 5351917a363SLisandro Dalcin /*@ 536bcf0153eSBarry Smith TSAdaptGetClip - Gets the admissible decrease/increase factor in step size in the time step adapter 5371917a363SLisandro Dalcin 5381917a363SLisandro Dalcin Not Collective 5391917a363SLisandro Dalcin 5404165533cSJose E. Roman Input Parameter: 5411917a363SLisandro Dalcin . adapt - adaptive controller context 5421917a363SLisandro Dalcin 5434165533cSJose E. Roman Output Parameters: 5441917a363SLisandro Dalcin + low - optional, admissible decrease factor 5451917a363SLisandro Dalcin - high - optional, admissible increase factor 5461917a363SLisandro Dalcin 5471917a363SLisandro Dalcin Level: intermediate 5481917a363SLisandro Dalcin 5491cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 5501917a363SLisandro Dalcin @*/ 551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high) 552d71ae5a4SJacob Faibussowitsch { 5531917a363SLisandro Dalcin PetscFunctionBegin; 5541917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 5554f572ea9SToby Isaac if (low) PetscAssertPointer(low, 2); 5564f572ea9SToby Isaac if (high) PetscAssertPointer(high, 3); 5571917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5581917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5601917a363SLisandro Dalcin } 5611917a363SLisandro Dalcin 5621917a363SLisandro Dalcin /*@ 563bcf0153eSBarry Smith TSAdaptSetScaleSolveFailed - Scale step size by this factor if solve fails 56462c23b28SMark 56520f4b53cSBarry Smith Logically Collective 56662c23b28SMark 5674165533cSJose E. Roman Input Parameters: 56862c23b28SMark + adapt - adaptive controller context 569ee300463SSatish Balay - scale - scale 57062c23b28SMark 571bcf0153eSBarry Smith Options Database Key: 57262c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 57362c23b28SMark 57462c23b28SMark Level: intermediate 57562c23b28SMark 5761cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 57762c23b28SMark @*/ 578d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale) 579d71ae5a4SJacob Faibussowitsch { 58062c23b28SMark PetscFunctionBegin; 58162c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 58262c23b28SMark PetscValidLogicalCollectiveReal(adapt, scale, 2); 58309cb0f53SBarry Smith PetscCheck(scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale); 58409cb0f53SBarry Smith PetscCheck(scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale); 58509cb0f53SBarry Smith adapt->scale_solve_failed = scale; 5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58762c23b28SMark } 58862c23b28SMark 58962c23b28SMark /*@ 59062c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 59162c23b28SMark 59262c23b28SMark Not Collective 59362c23b28SMark 5944165533cSJose E. Roman Input Parameter: 59562c23b28SMark . adapt - adaptive controller context 59662c23b28SMark 5974165533cSJose E. Roman Output Parameter: 598ee300463SSatish Balay . scale - scale factor 59962c23b28SMark 60062c23b28SMark Level: intermediate 60162c23b28SMark 6021cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 60362c23b28SMark @*/ 604d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale) 605d71ae5a4SJacob Faibussowitsch { 60662c23b28SMark PetscFunctionBegin; 60762c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 6084f572ea9SToby Isaac if (scale) PetscAssertPointer(scale, 2); 60962c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 6103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61162c23b28SMark } 61262c23b28SMark 61362c23b28SMark /*@ 614bcf0153eSBarry Smith TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the time step controller 6151917a363SLisandro Dalcin 61620f4b53cSBarry Smith Logically Collective 6171c3436cfSJed Brown 6184165533cSJose E. Roman Input Parameters: 619bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 6201c3436cfSJed Brown . hmin - minimum time step 6211c3436cfSJed Brown - hmax - maximum time step 6221c3436cfSJed Brown 6231c3436cfSJed Brown Options Database Keys: 624ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 625ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 6261c3436cfSJed Brown 6271c3436cfSJed Brown Level: intermediate 6281c3436cfSJed Brown 62909cb0f53SBarry Smith Note: 63009cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter 63109cb0f53SBarry Smith 63209cb0f53SBarry Smith Fortran Note: 63309cb0f53SBarry Smith Use `PETSC_CURRENT_REAL` 63409cb0f53SBarry Smith 6351cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 6361c3436cfSJed Brown @*/ 637d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax) 638d71ae5a4SJacob Faibussowitsch { 6391c3436cfSJed Brown PetscFunctionBegin; 6404782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 641b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmin, 2); 642b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmax, 3); 64309cb0f53SBarry 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); 64409cb0f53SBarry 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); 64509cb0f53SBarry Smith if (hmin != (PetscReal)PETSC_CURRENT) adapt->dt_min = hmin; 64609cb0f53SBarry Smith if (hmax != (PetscReal)PETSC_CURRENT) adapt->dt_max = hmax; 647b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 648b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6493c633725SBarry 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); 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6511917a363SLisandro Dalcin } 6521917a363SLisandro Dalcin 6531917a363SLisandro Dalcin /*@ 654bcf0153eSBarry Smith TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the time step controller 6551917a363SLisandro Dalcin 6561917a363SLisandro Dalcin Not Collective 6571917a363SLisandro Dalcin 6584165533cSJose E. Roman Input Parameter: 659bcf0153eSBarry Smith . adapt - time step adaptivity context, usually gotten with `TSGetAdapt()` 6601917a363SLisandro Dalcin 6614165533cSJose E. Roman Output Parameters: 6621917a363SLisandro Dalcin + hmin - minimum time step 6631917a363SLisandro Dalcin - hmax - maximum time step 6641917a363SLisandro Dalcin 6651917a363SLisandro Dalcin Level: intermediate 6661917a363SLisandro Dalcin 6671cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 6681917a363SLisandro Dalcin @*/ 669d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax) 670d71ae5a4SJacob Faibussowitsch { 6711917a363SLisandro Dalcin PetscFunctionBegin; 6721917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 6734f572ea9SToby Isaac if (hmin) PetscAssertPointer(hmin, 2); 6744f572ea9SToby Isaac if (hmax) PetscAssertPointer(hmax, 3); 6751917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6761917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6781c3436cfSJed Brown } 6791c3436cfSJed Brown 68014d0ab18SJacob Faibussowitsch /*@C 681bcf0153eSBarry Smith TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options. 68284df9cb4SJed Brown 683c3339decSBarry Smith Collective 68484df9cb4SJed Brown 68514d0ab18SJacob Faibussowitsch Input Parameters: 6862fe279fdSBarry Smith + adapt - the `TSAdapt` context 6872fe279fdSBarry Smith - PetscOptionsObject - object created by `PetscOptionsBegin()` 68884df9cb4SJed Brown 68984df9cb4SJed Brown Options Database Keys: 6901917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 691bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6921917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6931917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6941917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 69523de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 69623de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 69723de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 698de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 699de50f1caSBarry 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 70084df9cb4SJed Brown 70184df9cb4SJed Brown Level: advanced 70284df9cb4SJed Brown 703bcf0153eSBarry Smith Note: 704bcf0153eSBarry Smith This function is automatically called by `TSSetFromOptions()` 70584df9cb4SJed Brown 7061cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 707db781477SPatrick Sanan `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 70814d0ab18SJacob Faibussowitsch @*/ 709ce78bad3SBarry Smith PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems PetscOptionsObject) 710d71ae5a4SJacob Faibussowitsch { 71184df9cb4SJed Brown char type[256] = TSADAPTBASIC; 71262c23b28SMark PetscReal safety, reject_safety, clip[2], scale, hmin, hmax; 7131c3436cfSJed Brown PetscBool set, flg; 7141917a363SLisandro Dalcin PetscInt two; 71584df9cb4SJed Brown 71684df9cb4SJed Brown PetscFunctionBegin; 717dbbe0bcdSBarry Smith PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 71884df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 7191566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 720d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options"); 7219566063dSJacob 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)); 72248a46eb9SPierre Jolivet if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type)); 7231917a363SLisandro Dalcin 7249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set)); 7259566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg)); 7261917a363SLisandro Dalcin 7279371c9d4SSatish Balay safety = adapt->safety; 7289371c9d4SSatish Balay reject_safety = adapt->reject_safety; 7299566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set)); 7309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg)); 7319566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety)); 7321917a363SLisandro Dalcin 7339371c9d4SSatish Balay two = 2; 7349371c9d4SSatish Balay clip[0] = adapt->clip[0]; 7359371c9d4SSatish Balay clip[1] = adapt->clip[1]; 7369566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set)); 7373c633725SBarry Smith PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip"); 7389566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1])); 7391917a363SLisandro Dalcin 7409371c9d4SSatish Balay hmin = adapt->dt_min; 7419371c9d4SSatish Balay hmax = adapt->dt_max; 7429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set)); 7439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg)); 7449566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax)); 7451917a363SLisandro Dalcin 7469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set)); 7479566063dSJacob 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)); 748d580f011SEmil Constantinescu 7499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set)); 7509566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale)); 7511917a363SLisandro Dalcin 7529566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL)); 7533c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported"); 7541917a363SLisandro Dalcin 7559566063dSJacob 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)); 756de50f1caSBarry Smith 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set)); 7589566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt, flg)); 7591917a363SLisandro Dalcin 760dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 761d0609cedSBarry Smith PetscOptionsHeadEnd(); 7623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76384df9cb4SJed Brown } 76484df9cb4SJed Brown 76584df9cb4SJed Brown /*@ 76684df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 76784df9cb4SJed Brown 76820f4b53cSBarry Smith Logically Collective 76984df9cb4SJed Brown 7704165533cSJose E. Roman Input Parameter: 77184df9cb4SJed Brown . adapt - adaptive controller 77284df9cb4SJed Brown 77384df9cb4SJed Brown Level: developer 77484df9cb4SJed Brown 7751cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 77684df9cb4SJed Brown @*/ 777d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 778d71ae5a4SJacob Faibussowitsch { 77984df9cb4SJed Brown PetscFunctionBegin; 7804782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 7819566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates))); 7823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78384df9cb4SJed Brown } 78484df9cb4SJed Brown 78584df9cb4SJed Brown /*@C 78684df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 78784df9cb4SJed Brown 78820f4b53cSBarry Smith Logically Collective; No Fortran Support 78984df9cb4SJed Brown 7904165533cSJose E. Roman Input Parameters: 791bcf0153eSBarry Smith + adapt - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()` 79284df9cb4SJed Brown . name - name of the candidate scheme to add 79384df9cb4SJed Brown . order - order of the candidate scheme 79484df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7958d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 79684df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 79784df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 79884df9cb4SJed Brown 79984df9cb4SJed Brown Level: developer 80084df9cb4SJed Brown 8011cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 80284df9cb4SJed Brown @*/ 803d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse) 804d71ae5a4SJacob Faibussowitsch { 80584df9cb4SJed Brown PetscInt c; 80684df9cb4SJed Brown 80784df9cb4SJed Brown PetscFunctionBegin; 80884df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 80963a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order); 81084df9cb4SJed Brown if (inuse) { 8113c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 81284df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 81384df9cb4SJed Brown } 8141c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 8151c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 816bbd56ea5SKarl Rupp 81784df9cb4SJed Brown adapt->candidates.name[c] = name; 81884df9cb4SJed Brown adapt->candidates.order[c] = order; 81984df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 8208d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 82184df9cb4SJed Brown adapt->candidates.cost[c] = cost; 82284df9cb4SJed Brown adapt->candidates.n++; 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82484df9cb4SJed Brown } 82584df9cb4SJed Brown 8268d59e960SJed Brown /*@C 8278d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8288d59e960SJed Brown 8298d59e960SJed Brown Not Collective 8308d59e960SJed Brown 8314165533cSJose E. Roman Input Parameter: 8328d59e960SJed Brown . adapt - time step adaptivity context 8338d59e960SJed Brown 8344165533cSJose E. Roman Output Parameters: 8358d59e960SJed Brown + n - number of candidate schemes, always at least 1 8368d59e960SJed Brown . order - the order of each candidate scheme 8378d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8388d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8398d59e960SJed Brown - cost - the relative cost of each scheme 8408d59e960SJed Brown 8418d59e960SJed Brown Level: developer 8428d59e960SJed Brown 8438d59e960SJed Brown Note: 8448d59e960SJed Brown The current scheme is always returned in the first slot 8458d59e960SJed Brown 8461cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 8478d59e960SJed Brown @*/ 848d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost) 849d71ae5a4SJacob Faibussowitsch { 8508d59e960SJed Brown PetscFunctionBegin; 8518d59e960SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 8528d59e960SJed Brown if (n) *n = adapt->candidates.n; 8538d59e960SJed Brown if (order) *order = adapt->candidates.order; 8548d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8558d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8568d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8588d59e960SJed Brown } 8598d59e960SJed Brown 86084df9cb4SJed Brown /*@C 86184df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 86284df9cb4SJed Brown 863c3339decSBarry Smith Collective 86484df9cb4SJed Brown 8654165533cSJose E. Roman Input Parameters: 866da81f932SPierre Jolivet + adapt - adaptive controller 86797bb3fdcSJose E. Roman . ts - time stepper 86884df9cb4SJed Brown - h - current step size 86984df9cb4SJed Brown 8704165533cSJose E. Roman Output Parameters: 8711566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 87284df9cb4SJed Brown . next_h - step size to use for the next step 873bcf0153eSBarry Smith - accept - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size 8741c3436cfSJed Brown 87584df9cb4SJed Brown Level: developer 87684df9cb4SJed Brown 877bcf0153eSBarry Smith Note: 878bcf0153eSBarry Smith The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is 879bcf0153eSBarry Smith being retried after an initial rejection. 880bcf0153eSBarry Smith 8811cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 88284df9cb4SJed Brown @*/ 883d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept) 884d71ae5a4SJacob Faibussowitsch { 8851566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8861566a47fSLisandro Dalcin PetscInt scheme = 0; 8870b99f514SJed Brown PetscReal wlte = -1.0; 8885e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8895e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 89084df9cb4SJed Brown 89184df9cb4SJed Brown PetscFunctionBegin; 89284df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 89384df9cb4SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 8944f572ea9SToby Isaac if (next_sc) PetscAssertPointer(next_sc, 4); 8954f572ea9SToby Isaac PetscAssertPointer(next_h, 5); 8964f572ea9SToby Isaac PetscAssertPointer(accept, 6); 8971566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8981566a47fSLisandro Dalcin 8991566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events */ 900ca4445c7SIlya Fursov if (ts->event && ts->event->processing) { 9011566a47fSLisandro Dalcin *next_h = h; 9021566a47fSLisandro Dalcin *accept = PETSC_TRUE; 903fe4ad979SIlya Fursov if (adapt->monitor) { 904fe4ad979SIlya Fursov PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 905fe4ad979SIlya Fursov 906fe4ad979SIlya Fursov if (ts->event->iterctr == 0) { 907fe4ad979SIlya Fursov /* 908fe4ad979SIlya Fursov An event has been found, now finalising the event processing: performing the 1st and 2nd post-event steps. 909fe4ad979SIlya Fursov Entering this if-branch means both these steps (set to either PETSC_DECIDE or numerical value) are managed 910fe4ad979SIlya Fursov by the event handler. In this case the 1st post-event step is always accepted, without interference of TSAdapt. 911fe4ad979SIlya Fursov Note: if the 2nd post-event step is not managed by the event handler (e.g. given 1st = numerical, 2nd = PETSC_DECIDE), 912fe4ad979SIlya Fursov this if-branch is not entered, and TSAdapt may reject/adjust the proposed 1st post-event step. 913fe4ad979SIlya Fursov */ 914fe4ad979SIlya 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)); 915fe4ad979SIlya Fursov } else PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt does not interfere, step %3" PetscInt_FMT " accepted. Event handling in progress\n", ts->steps)); 916fe4ad979SIlya Fursov 917fe4ad979SIlya Fursov PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 918fe4ad979SIlya Fursov } 9193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9201566a47fSLisandro Dalcin } 9211566a47fSLisandro Dalcin 922dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter); 92363a3b9bcSJacob 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); 9243c633725SBarry Smith PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h); 9251566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 9261566a47fSLisandro Dalcin 92749354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 92836b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 929ca4445c7SIlya Fursov PetscReal t = ts->ptime + ts->time_step, tend, tmax, h1, hmax; 930fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 931fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 932ca4445c7SIlya Fursov PetscReal eps = 10 * PETSC_MACHINE_EPSILON; 9334a658b32SHong Zhang 934fe4ad979SIlya Fursov /* 935fe4ad979SIlya Fursov Logic in using 'dt_span_cached': 936fe4ad979SIlya Fursov 1. It always overrides *next_h, except (any of): 937fe4ad979SIlya Fursov a) the current step was rejected, 938fe4ad979SIlya Fursov b) the adaptor proposed to decrease the next step, 939fe4ad979SIlya Fursov c) the adaptor proposed *next_h > dt_span_cached. 940136cf249SJames Wright 2. If *next_h was adjusted by eval_times points (or the final point): 941fe4ad979SIlya Fursov -- when dt_span_cached is filled (>0), it keeps its value, 942fe4ad979SIlya Fursov -- when dt_span_cached is clear (==0), it gets the unadjusted version of *next_h. 943fe4ad979SIlya Fursov 3. If *next_h was not adjusted as in (2), dt_span_cached is cleared. 944fe4ad979SIlya Fursov Note, if a combination (1.b || 1.c) && (3) takes place, this means that 945fe4ad979SIlya Fursov dt_span_cached remains unused at the moment of clearing. 946fe4ad979SIlya Fursov If (1.a) takes place, dt_span_cached keeps its value. 947fe4ad979SIlya Fursov Also, dt_span_cached can be updated by the event handler, see tsevent.c. 948fe4ad979SIlya Fursov */ 949136cf249SJames Wright if (h <= *next_h && *next_h <= adapt->dt_eval_times_cached) *next_h = adapt->dt_eval_times_cached; /* try employing the cache */ 950ca4445c7SIlya Fursov h1 = *next_h; 951ca4445c7SIlya Fursov tend = t + h1; 952ca4445c7SIlya Fursov 953136cf249SJames Wright if (ts->eval_times && ts->eval_times->time_point_idx < ts->eval_times->num_time_points) { 954136cf249SJames Wright PetscCheck(ts->eval_times->worktol == 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_PLIB, "Unexpected state (tspan->worktol != 0) in TSAdaptChoose()"); 955136cf249SJames Wright ts->eval_times->worktol = ts->eval_times->reltol * h1 + ts->eval_times->abstol; 956136cf249SJames 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 */ 957136cf249SJames 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]; 9584a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */ 959136cf249SJames Wright else tmax = ts->eval_times->time_points[ts->eval_times->time_point_idx]; 9604a658b32SHong Zhang } else tmax = ts->max_time; 961ca4445c7SIlya Fursov tmax = PetscMin(tmax, ts->max_time); 9624a658b32SHong Zhang hmax = tmax - t; 963ca4445c7SIlya 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()"); 964ca4445c7SIlya Fursov 96536b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 966ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * b > hmax) *next_h = hmax / 2; 967ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * a > hmax) *next_h = hmax; 968136cf249SJames 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 */ 969136cf249SJames 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 */ 970e1db57b0SHong Zhang } 9711c3436cfSJed Brown if (adapt->monitor) { 9721566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 9739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9740b99f514SJed Brown if (wlte < 0) { 9759371c9d4SSatish 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", 9769371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h)); 9770b99f514SJed Brown } else { 9789371c9d4SSatish 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", 9799371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter)); 9800b99f514SJed Brown } 9819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9821c3436cfSJed Brown } 9833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98484df9cb4SJed Brown } 98584df9cb4SJed Brown 98697335746SJed Brown /*@ 987de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 988de50f1caSBarry Smith before increasing the time step. 989de50f1caSBarry Smith 990c3339decSBarry Smith Logicially Collective 991de50f1caSBarry Smith 9924165533cSJose E. Roman Input Parameters: 993de50f1caSBarry Smith + adapt - adaptive controller context 994de50f1caSBarry Smith - cnt - the number of timesteps 995de50f1caSBarry Smith 996de50f1caSBarry Smith Options Database Key: 997de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 998de50f1caSBarry Smith 999de50f1caSBarry Smith Level: advanced 1000de50f1caSBarry Smith 1001bcf0153eSBarry Smith Notes: 1002bcf0153eSBarry Smith This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 1003bcf0153eSBarry Smith 1004bcf0153eSBarry Smith The successful use of this option is problem dependent 1005bcf0153eSBarry Smith 1006b43aa488SJacob Faibussowitsch Developer Notes: 1007bcf0153eSBarry Smith There is no theory to support this option 1008bcf0153eSBarry Smith 10091cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt` 1010de50f1caSBarry Smith @*/ 1011d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt) 1012d71ae5a4SJacob Faibussowitsch { 1013de50f1caSBarry Smith PetscFunctionBegin; 1014de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 10153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1016de50f1caSBarry Smith } 1017de50f1caSBarry Smith 1018de50f1caSBarry Smith /*@ 10196bc98fa9SBarry 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) 102097335746SJed Brown 1021c3339decSBarry Smith Collective 102297335746SJed Brown 10234165533cSJose E. Roman Input Parameters: 102497335746SJed Brown + adapt - adaptive controller context 1025b295832fSPierre Barbier de Reuille . ts - time stepper 1026b295832fSPierre Barbier de Reuille . t - Current simulation time 1027b295832fSPierre Barbier de Reuille - Y - Current solution vector 102897335746SJed Brown 10294165533cSJose E. Roman Output Parameter: 1030bcf0153eSBarry Smith . accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject 103197335746SJed Brown 103297335746SJed Brown Level: developer 103397335746SJed Brown 10341cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt` 103597335746SJed Brown @*/ 1036d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) 1037d71ae5a4SJacob Faibussowitsch { 10381566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 10396c6709e3SStefano Zampini PetscBool func_accept, snes_div_func; 104097335746SJed Brown 104197335746SJed Brown PetscFunctionBegin; 10424782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 10434782b174SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 10444f572ea9SToby Isaac PetscAssertPointer(accept, 5); 10451566a47fSLisandro Dalcin 10466c6709e3SStefano Zampini PetscCall(TSFunctionDomainError(ts, t, Y, &func_accept)); 10479566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 10486c6709e3SStefano Zampini snes_div_func = (PetscBool)(snesreason == SNES_DIVERGED_FUNCTION_DOMAIN); 10496c6709e3SStefano Zampini if (func_accept && snesreason < 0 && !snes_div_func) { 105097335746SJed Brown *accept = PETSC_FALSE; 10516c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failure: %s\n", ts->steps, SNESConvergedReasons[snesreason])); 105209cb0f53SBarry Smith if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures != PETSC_UNLIMITED) { 105397335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 105463a3b9bcSJacob 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)); 105597335746SJed Brown if (adapt->monitor) { 10569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10579371c9d4SSatish 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, 10589371c9d4SSatish Balay (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures)); 10599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 106097335746SJed Brown } 1061cb9d8021SPierre Barbier de Reuille } 1062cb9d8021SPierre Barbier de Reuille } else { 10636c6709e3SStefano Zampini *accept = (PetscBool)(func_accept && !snes_div_func); 10646c6709e3SStefano Zampini if (*accept && adapt->checkstage) PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept)); 10656bc98fa9SBarry Smith if (!*accept) { 10666c6709e3SStefano Zampini const char *user_func = !func_accept ? "TSSetFunctionDomainError()" : "TSAdaptSetCheckStage"; 10676c6709e3SStefano Zampini const char *snes_err = "SNES invalid function domain"; 10686c6709e3SStefano Zampini const char *err_msg = snes_div_func && func_accept ? snes_err : user_func; 10696c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by %s\n", ts->steps, err_msg)); 10706bc98fa9SBarry Smith if (adapt->monitor) { 10719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10726c6709e3SStefano Zampini PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected by %s\n", ((PetscObject)adapt)->type_name, ts->steps, err_msg)); 10739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10746bc98fa9SBarry Smith } 10756bc98fa9SBarry Smith } 1076cb9d8021SPierre Barbier de Reuille } 1077cb9d8021SPierre Barbier de Reuille 107857508eceSPierre Jolivet if (!*accept && !ts->reason) { 10791566a47fSLisandro Dalcin PetscReal dt, new_dt; 10809566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1081cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 10829566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, new_dt)); 1083de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 108497335746SJed Brown if (adapt->monitor) { 10859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 10866c6709e3SStefano 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], 10876c6709e3SStefano Zampini (double)ts->ptime, (double)dt, (double)new_dt)); 10889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 108997335746SJed Brown } 109097335746SJed Brown } 10913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109297335746SJed Brown } 109397335746SJed Brown 109484df9cb4SJed Brown /*@ 109584df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 109684df9cb4SJed Brown 1097d083f849SBarry Smith Collective 109884df9cb4SJed Brown 109984df9cb4SJed Brown Input Parameter: 110084df9cb4SJed Brown . comm - The communicator 110184df9cb4SJed Brown 110284df9cb4SJed Brown Output Parameter: 1103b43aa488SJacob Faibussowitsch . inadapt - new `TSAdapt` object 110484df9cb4SJed Brown 110584df9cb4SJed Brown Level: developer 110684df9cb4SJed Brown 1107bcf0153eSBarry Smith Note: 1108bcf0153eSBarry Smith `TSAdapt` creation is handled by `TS`, so users should not need to call this function. 110984df9cb4SJed Brown 11101cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 111184df9cb4SJed Brown @*/ 1112d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt) 1113d71ae5a4SJacob Faibussowitsch { 111484df9cb4SJed Brown TSAdapt adapt; 111584df9cb4SJed Brown 111684df9cb4SJed Brown PetscFunctionBegin; 11174f572ea9SToby Isaac PetscAssertPointer(inadapt, 2); 11189566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 11193b3bcf4cSLisandro Dalcin 11209566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView)); 1121bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 11221917a363SLisandro Dalcin adapt->safety = 0.9; 11231917a363SLisandro Dalcin adapt->reject_safety = 0.5; 11241917a363SLisandro Dalcin adapt->clip[0] = 0.1; 11251917a363SLisandro Dalcin adapt->clip[1] = 10.; 11261c3436cfSJed Brown adapt->dt_min = 1e-20; 11271917a363SLisandro Dalcin adapt->dt_max = 1e+20; 11281c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1129d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 113097335746SJed Brown adapt->scale_solve_failed = 0.25; 1131fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1132fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1133fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1134fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 11357619abb3SShri adapt->wnormtype = NORM_2; 1136de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 113784df9cb4SJed Brown *inadapt = adapt; 11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113984df9cb4SJed Brown } 1140