184df9cb4SJed Brown 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 384df9cb4SJed Brown 41b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID; 51b9b13dfSLisandro Dalcin 6140e18c1SBarry Smith static PetscFunctionList TSAdaptList; 784df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 884df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 984df9cb4SJed Brown 108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 111566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 12cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt); 138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 141917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 15d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt); 1684df9cb4SJed Brown 1784df9cb4SJed Brown /*@C 181c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation 191c84c290SBarry Smith 201c84c290SBarry Smith Not Collective 211c84c290SBarry Smith 221c84c290SBarry Smith Input Parameters: 231c84c290SBarry Smith + name_scheme - name of user-defined adaptivity scheme 241c84c290SBarry Smith - routine_create - routine to create method context 251c84c290SBarry Smith 261c84c290SBarry Smith Notes: 271c84c290SBarry Smith TSAdaptRegister() may be called multiple times to add several user-defined families. 281c84c290SBarry Smith 291c84c290SBarry Smith Sample usage: 301c84c290SBarry Smith .vb 31bdf89e91SBarry Smith TSAdaptRegister("my_scheme",MySchemeCreate); 321c84c290SBarry Smith .ve 331c84c290SBarry Smith 341c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via 351c84c290SBarry Smith $ TSAdaptSetType(ts,"my_scheme") 361c84c290SBarry Smith or at runtime via the option 371c84c290SBarry Smith $ -ts_adapt_type my_scheme 3884df9cb4SJed Brown 3984df9cb4SJed Brown Level: advanced 401c84c290SBarry Smith 41db781477SPatrick Sanan .seealso: `TSAdaptRegisterAll()` 4284df9cb4SJed Brown @*/ 439371c9d4SSatish Balay PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt)) { 4484df9cb4SJed Brown PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 469566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function)); 4784df9cb4SJed Brown PetscFunctionReturn(0); 4884df9cb4SJed Brown } 4984df9cb4SJed Brown 5084df9cb4SJed Brown /*@C 5184df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5284df9cb4SJed Brown 5384df9cb4SJed Brown Not Collective 5484df9cb4SJed Brown 5584df9cb4SJed Brown Level: advanced 5684df9cb4SJed Brown 57db781477SPatrick Sanan .seealso: `TSAdaptRegisterDestroy()` 5884df9cb4SJed Brown @*/ 599371c9d4SSatish Balay PetscErrorCode TSAdaptRegisterAll(void) { 6084df9cb4SJed Brown PetscFunctionBegin; 610f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 620f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 639566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None)); 649566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic)); 659566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP)); 669566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL)); 679566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE)); 689566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History)); 6984df9cb4SJed Brown PetscFunctionReturn(0); 7084df9cb4SJed Brown } 7184df9cb4SJed Brown 7284df9cb4SJed Brown /*@C 73560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 7484df9cb4SJed Brown called from PetscFinalize(). 7584df9cb4SJed Brown 7684df9cb4SJed Brown Level: developer 7784df9cb4SJed Brown 78db781477SPatrick Sanan .seealso: `PetscFinalize()` 7984df9cb4SJed Brown @*/ 809371c9d4SSatish Balay PetscErrorCode TSAdaptFinalizePackage(void) { 8184df9cb4SJed Brown PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSAdaptList)); 8384df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 8484df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 8584df9cb4SJed Brown PetscFunctionReturn(0); 8684df9cb4SJed Brown } 8784df9cb4SJed Brown 8884df9cb4SJed Brown /*@C 8984df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 908a690491SBarry Smith called from TSInitializePackage(). 9184df9cb4SJed Brown 9284df9cb4SJed Brown Level: developer 9384df9cb4SJed Brown 94db781477SPatrick Sanan .seealso: `PetscInitialize()` 9584df9cb4SJed Brown @*/ 969371c9d4SSatish Balay PetscErrorCode TSAdaptInitializePackage(void) { 9784df9cb4SJed Brown PetscFunctionBegin; 9884df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 9984df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 1009566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID)); 1019566063dSJacob Faibussowitsch PetscCall(TSAdaptRegisterAll()); 1029566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage)); 10384df9cb4SJed Brown PetscFunctionReturn(0); 10484df9cb4SJed Brown } 10584df9cb4SJed Brown 1067eef6055SBarry Smith /*@C 1077eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1087eef6055SBarry Smith 1097eef6055SBarry Smith Logicially Collective on TSAdapt 1107eef6055SBarry Smith 111d8d19677SJose E. Roman Input Parameters: 112d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1137eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1147eef6055SBarry Smith 1157eef6055SBarry Smith Options Database: 116ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1177eef6055SBarry Smith 1187eef6055SBarry Smith Level: intermediate 1197eef6055SBarry Smith 120db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()` 1217eef6055SBarry Smith @*/ 1229371c9d4SSatish Balay PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type) { 123ef749922SLisandro Dalcin PetscBool match; 1245f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSAdapt); 12584df9cb4SJed Brown 12684df9cb4SJed Brown PetscFunctionBegin; 1274782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 128b92453a8SLisandro Dalcin PetscValidCharPointer(type, 2); 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match)); 130ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSAdaptList, type, &r)); 1323c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type); 133dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, destroy); 1349566063dSJacob Faibussowitsch PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps))); 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 1369566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 13784df9cb4SJed Brown PetscFunctionReturn(0); 13884df9cb4SJed Brown } 13984df9cb4SJed Brown 140d0288e4fSLisandro Dalcin /*@C 141d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 142d0288e4fSLisandro Dalcin 143d0288e4fSLisandro Dalcin Not Collective 144d0288e4fSLisandro Dalcin 145d0288e4fSLisandro Dalcin Input Parameter: 146d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 147d0288e4fSLisandro Dalcin 148d0288e4fSLisandro Dalcin Output Parameter: 149d0288e4fSLisandro Dalcin . type - The name of TS adapter method 150d0288e4fSLisandro Dalcin 151d0288e4fSLisandro Dalcin Level: intermediate 152d0288e4fSLisandro Dalcin 153db781477SPatrick Sanan .seealso `TSAdaptSetType()` 154d0288e4fSLisandro Dalcin @*/ 1559371c9d4SSatish Balay PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type) { 156d0288e4fSLisandro Dalcin PetscFunctionBegin; 157d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 158d0288e4fSLisandro Dalcin PetscValidPointer(type, 2); 159d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 160d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 161d0288e4fSLisandro Dalcin } 162d0288e4fSLisandro Dalcin 1639371c9d4SSatish Balay PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[]) { 16484df9cb4SJed Brown PetscFunctionBegin; 1654782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 16784df9cb4SJed Brown PetscFunctionReturn(0); 16884df9cb4SJed Brown } 16984df9cb4SJed Brown 170ad6bc421SBarry Smith /*@C 171ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 172ad6bc421SBarry Smith 173ad6bc421SBarry Smith Collective on PetscViewer 174ad6bc421SBarry Smith 175ad6bc421SBarry Smith Input Parameters: 176ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 177ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 178ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 179ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 180ad6bc421SBarry Smith 181ad6bc421SBarry Smith Level: intermediate 182ad6bc421SBarry Smith 183ad6bc421SBarry Smith Notes: 184ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 185ad6bc421SBarry Smith 186ad6bc421SBarry Smith Notes for advanced users: 187ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 188ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 189ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 190ad6bc421SBarry Smith format is 191ad6bc421SBarry Smith .vb 192ad6bc421SBarry Smith has not yet been determined 193ad6bc421SBarry Smith .ve 194ad6bc421SBarry Smith 195db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()` 196ad6bc421SBarry Smith @*/ 1979371c9d4SSatish Balay PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer) { 198ad6bc421SBarry Smith PetscBool isbinary; 199ad6bc421SBarry Smith char type[256]; 200ad6bc421SBarry Smith 201ad6bc421SBarry Smith PetscFunctionBegin; 2024782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 203ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 2053c633725SBarry Smith PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 206ad6bc421SBarry Smith 2079566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 2089566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, type)); 209dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, load, viewer); 210ad6bc421SBarry Smith PetscFunctionReturn(0); 211ad6bc421SBarry Smith } 212ad6bc421SBarry Smith 2139371c9d4SSatish Balay PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer) { 2141c167fc2SEmil Constantinescu PetscBool iascii, isbinary, isnone, isglee; 21584df9cb4SJed Brown 21684df9cb4SJed Brown PetscFunctionBegin; 2174782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 2189566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer)); 2194782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2204782b174SLisandro Dalcin PetscCheckSameComm(adapt, 1, viewer, 2); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 22384df9cb4SJed Brown if (iascii) { 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee)); 2271917a363SLisandro Dalcin if (!isnone) { 2289566063dSJacob Faibussowitsch if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, " always accepting steps\n")); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " safety factor %g\n", (double)adapt->safety)); 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " extra safety factor after step rejection %g\n", (double)adapt->reject_safety)); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest increase %g\n", (double)adapt->clip[1])); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest decrease %g\n", (double)adapt->clip[0])); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum allowed timestep %g\n", (double)adapt->dt_max)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " minimum allowed timestep %g\n", (double)adapt->dt_min)); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max)); 2361c167fc2SEmil Constantinescu } 2371c167fc2SEmil Constantinescu if (isglee) { 2381c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses local error control\n")); 2401c167fc2SEmil Constantinescu } else { 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses global error control\n")); 2421c167fc2SEmil Constantinescu } 2431917a363SLisandro Dalcin } 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 245dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, view, viewer); 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 247ad6bc421SBarry Smith } else if (isbinary) { 248ad6bc421SBarry Smith char type[256]; 249ad6bc421SBarry Smith 250ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 2519566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 253dbbe0bcdSBarry Smith } else PetscTryTypeMethod(adapt, view, viewer); 25484df9cb4SJed Brown PetscFunctionReturn(0); 25584df9cb4SJed Brown } 25684df9cb4SJed Brown 257099af0a3SLisandro Dalcin /*@ 258099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 259099af0a3SLisandro Dalcin 260099af0a3SLisandro Dalcin Collective on TS 261099af0a3SLisandro Dalcin 262099af0a3SLisandro Dalcin Input Parameter: 263099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 264099af0a3SLisandro Dalcin 265099af0a3SLisandro Dalcin Level: developer 266099af0a3SLisandro Dalcin 267db781477SPatrick Sanan .seealso: `TSAdaptCreate()`, `TSAdaptDestroy()` 268099af0a3SLisandro Dalcin @*/ 2699371c9d4SSatish Balay PetscErrorCode TSAdaptReset(TSAdapt adapt) { 270099af0a3SLisandro Dalcin PetscFunctionBegin; 271099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 272dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, reset); 273099af0a3SLisandro Dalcin PetscFunctionReturn(0); 274099af0a3SLisandro Dalcin } 275099af0a3SLisandro Dalcin 2769371c9d4SSatish Balay PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) { 27784df9cb4SJed Brown PetscFunctionBegin; 27884df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 27984df9cb4SJed Brown PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1); 2809371c9d4SSatish Balay if (--((PetscObject)(*adapt))->refct > 0) { 2819371c9d4SSatish Balay *adapt = NULL; 2829371c9d4SSatish Balay PetscFunctionReturn(0); 2839371c9d4SSatish Balay } 284099af0a3SLisandro Dalcin 2859566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt)); 286099af0a3SLisandro Dalcin 287dbbe0bcdSBarry Smith PetscTryTypeMethod((*adapt), destroy); 2889566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 2899566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 29084df9cb4SJed Brown PetscFunctionReturn(0); 29184df9cb4SJed Brown } 29284df9cb4SJed Brown 2931c3436cfSJed Brown /*@ 2941c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2951c3436cfSJed Brown 2961c3436cfSJed Brown Collective on TSAdapt 2971c3436cfSJed Brown 2984165533cSJose E. Roman Input Parameters: 2991c3436cfSJed Brown + adapt - adaptive controller context 3001c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3011c3436cfSJed Brown 302bf997491SLisandro Dalcin Options Database Keys: 303ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 304bf997491SLisandro Dalcin 3051c3436cfSJed Brown Level: intermediate 3061c3436cfSJed Brown 307db781477SPatrick Sanan .seealso: `TSAdaptChoose()` 3081c3436cfSJed Brown @*/ 3099371c9d4SSatish Balay PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg) { 3101c3436cfSJed Brown PetscFunctionBegin; 3114782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3124782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flg, 2); 3131c3436cfSJed Brown if (flg) { 3149566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor)); 3151c3436cfSJed Brown } else { 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor)); 3171c3436cfSJed Brown } 3181c3436cfSJed Brown PetscFunctionReturn(0); 3191c3436cfSJed Brown } 3201c3436cfSJed Brown 3210873808bSJed Brown /*@C 322bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3230873808bSJed Brown 3240873808bSJed Brown Logically collective on TSAdapt 3250873808bSJed Brown 3264165533cSJose E. Roman Input Parameters: 3270873808bSJed Brown + adapt - adaptive controller context 3280873808bSJed Brown - func - stage check function 3290873808bSJed Brown 3300873808bSJed Brown Arguments of func: 3310873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3320873808bSJed Brown 3330873808bSJed Brown + adapt - adaptive controller context 3340873808bSJed Brown . ts - time stepping context 3350873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3360873808bSJed Brown 3370873808bSJed Brown Level: advanced 3380873808bSJed Brown 339db781477SPatrick Sanan .seealso: `TSAdaptChoose()` 3400873808bSJed Brown @*/ 3419371c9d4SSatish Balay PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt, TS, PetscReal, Vec, PetscBool *)) { 3420873808bSJed Brown PetscFunctionBegin; 3430873808bSJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 34468ae941aSLisandro Dalcin adapt->checkstage = func; 3450873808bSJed Brown PetscFunctionReturn(0); 3460873808bSJed Brown } 3470873808bSJed Brown 3481c3436cfSJed Brown /*@ 349bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 350bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 351bf997491SLisandro Dalcin 3521917a363SLisandro Dalcin Logically collective on TSAdapt 353bf997491SLisandro Dalcin 3544165533cSJose E. Roman Input Parameters: 355bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 356bf997491SLisandro Dalcin - flag - whether to always accept steps 357bf997491SLisandro Dalcin 358bf997491SLisandro Dalcin Options Database Keys: 359ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 360bf997491SLisandro Dalcin 361bf997491SLisandro Dalcin Level: intermediate 362bf997491SLisandro Dalcin 363db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptChoose()` 364bf997491SLisandro Dalcin @*/ 3659371c9d4SSatish Balay PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag) { 366bf997491SLisandro Dalcin PetscFunctionBegin; 367bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 368bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flag, 2); 369bf997491SLisandro Dalcin adapt->always_accept = flag; 370bf997491SLisandro Dalcin PetscFunctionReturn(0); 371bf997491SLisandro Dalcin } 372bf997491SLisandro Dalcin 373bf997491SLisandro Dalcin /*@ 3741917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 3751c3436cfSJed Brown 3761917a363SLisandro Dalcin Logically collective on TSAdapt 3771917a363SLisandro Dalcin 3784165533cSJose E. Roman Input Parameters: 3791917a363SLisandro Dalcin + adapt - adaptive controller context 3801917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 381ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 3821917a363SLisandro Dalcin 3831917a363SLisandro Dalcin Options Database Keys: 384ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 385ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 3861917a363SLisandro Dalcin 3871917a363SLisandro Dalcin Level: intermediate 3881917a363SLisandro Dalcin 389db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 3901917a363SLisandro Dalcin @*/ 3919371c9d4SSatish Balay PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety) { 3921917a363SLisandro Dalcin PetscFunctionBegin; 3931917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3941917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, safety, 2); 3951917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, reject_safety, 3); 3963c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety); 3973c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety); 3983c633725SBarry Smith PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be non negative", (double)reject_safety); 3993c633725SBarry Smith PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be less than one", (double)reject_safety); 4001917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4011917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4021917a363SLisandro Dalcin PetscFunctionReturn(0); 4031917a363SLisandro Dalcin } 4041917a363SLisandro Dalcin 4051917a363SLisandro Dalcin /*@ 4061917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4071917a363SLisandro Dalcin 4081917a363SLisandro Dalcin Not Collective 4091917a363SLisandro Dalcin 4104165533cSJose E. Roman Input Parameter: 4111917a363SLisandro Dalcin . adapt - adaptive controller context 4121917a363SLisandro Dalcin 4134165533cSJose E. Roman Output Parameters: 4141917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4151917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4161917a363SLisandro Dalcin 4171917a363SLisandro Dalcin Level: intermediate 4181917a363SLisandro Dalcin 419db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 4201917a363SLisandro Dalcin @*/ 4219371c9d4SSatish Balay PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety) { 4221917a363SLisandro Dalcin PetscFunctionBegin; 4231917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4241917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety, 2); 4251917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety, 3); 4261917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4271917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4281917a363SLisandro Dalcin PetscFunctionReturn(0); 4291917a363SLisandro Dalcin } 4301917a363SLisandro Dalcin 4311917a363SLisandro Dalcin /*@ 43276cddca1SEmil Constantinescu TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components. 43376cddca1SEmil Constantinescu 43476cddca1SEmil Constantinescu Logically collective on TSAdapt 43576cddca1SEmil Constantinescu 4364165533cSJose E. Roman Input Parameters: 43776cddca1SEmil Constantinescu + adapt - adaptive controller context 43876cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 43976cddca1SEmil Constantinescu 44076cddca1SEmil Constantinescu Options Database Keys: 44176cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 44276cddca1SEmil Constantinescu 44376cddca1SEmil Constantinescu Level: intermediate 44476cddca1SEmil Constantinescu 445db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 44676cddca1SEmil Constantinescu @*/ 4479371c9d4SSatish Balay PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore) { 44876cddca1SEmil Constantinescu PetscFunctionBegin; 44976cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 45076cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt, max_ignore, 2); 45176cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 45276cddca1SEmil Constantinescu PetscFunctionReturn(0); 45376cddca1SEmil Constantinescu } 45476cddca1SEmil Constantinescu 45576cddca1SEmil Constantinescu /*@ 45676cddca1SEmil Constantinescu TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). 45776cddca1SEmil Constantinescu 45876cddca1SEmil Constantinescu Not Collective 45976cddca1SEmil Constantinescu 4604165533cSJose E. Roman Input Parameter: 46176cddca1SEmil Constantinescu . adapt - adaptive controller context 46276cddca1SEmil Constantinescu 4634165533cSJose E. Roman Output Parameter: 46476cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 46576cddca1SEmil Constantinescu 46676cddca1SEmil Constantinescu Level: intermediate 46776cddca1SEmil Constantinescu 468db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 46976cddca1SEmil Constantinescu @*/ 4709371c9d4SSatish Balay PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore) { 47176cddca1SEmil Constantinescu PetscFunctionBegin; 47276cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 47376cddca1SEmil Constantinescu PetscValidRealPointer(max_ignore, 2); 47476cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 47576cddca1SEmil Constantinescu PetscFunctionReturn(0); 47676cddca1SEmil Constantinescu } 47776cddca1SEmil Constantinescu 47876cddca1SEmil Constantinescu /*@ 4791917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4801917a363SLisandro Dalcin 4811917a363SLisandro Dalcin Logically collective on TSAdapt 4821917a363SLisandro Dalcin 4834165533cSJose E. Roman Input Parameters: 4841917a363SLisandro Dalcin + adapt - adaptive controller context 4851917a363SLisandro Dalcin . low - admissible decrease factor 486ec18b7bcSLisandro Dalcin - high - admissible increase factor 4871917a363SLisandro Dalcin 4881917a363SLisandro Dalcin Options Database Keys: 489ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 4901917a363SLisandro Dalcin 4911917a363SLisandro Dalcin Level: intermediate 4921917a363SLisandro Dalcin 493db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 4941917a363SLisandro Dalcin @*/ 4959371c9d4SSatish Balay PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high) { 4961917a363SLisandro Dalcin PetscFunctionBegin; 4971917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 4981917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, low, 2); 4991917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, high, 3); 5003c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low); 5013c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low); 5023c633725SBarry Smith PetscCheck(high == PETSC_DEFAULT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high); 5031917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 5041917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 5051917a363SLisandro Dalcin PetscFunctionReturn(0); 5061917a363SLisandro Dalcin } 5071917a363SLisandro Dalcin 5081917a363SLisandro Dalcin /*@ 5091917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5101917a363SLisandro Dalcin 5111917a363SLisandro Dalcin Not Collective 5121917a363SLisandro Dalcin 5134165533cSJose E. Roman Input Parameter: 5141917a363SLisandro Dalcin . adapt - adaptive controller context 5151917a363SLisandro Dalcin 5164165533cSJose E. Roman Output Parameters: 5171917a363SLisandro Dalcin + low - optional, admissible decrease factor 5181917a363SLisandro Dalcin - high - optional, admissible increase factor 5191917a363SLisandro Dalcin 5201917a363SLisandro Dalcin Level: intermediate 5211917a363SLisandro Dalcin 522db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 5231917a363SLisandro Dalcin @*/ 5249371c9d4SSatish Balay PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high) { 5251917a363SLisandro Dalcin PetscFunctionBegin; 5261917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 5271917a363SLisandro Dalcin if (low) PetscValidRealPointer(low, 2); 5281917a363SLisandro Dalcin if (high) PetscValidRealPointer(high, 3); 5291917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5301917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5311917a363SLisandro Dalcin PetscFunctionReturn(0); 5321917a363SLisandro Dalcin } 5331917a363SLisandro Dalcin 5341917a363SLisandro Dalcin /*@ 53562c23b28SMark TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 53662c23b28SMark 53762c23b28SMark Logically collective on TSAdapt 53862c23b28SMark 5394165533cSJose E. Roman Input Parameters: 54062c23b28SMark + adapt - adaptive controller context 541ee300463SSatish Balay - scale - scale 54262c23b28SMark 54362c23b28SMark Options Database Keys: 54462c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 54562c23b28SMark 54662c23b28SMark Level: intermediate 54762c23b28SMark 548db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 54962c23b28SMark @*/ 5509371c9d4SSatish Balay PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale) { 55162c23b28SMark PetscFunctionBegin; 55262c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 55362c23b28SMark PetscValidLogicalCollectiveReal(adapt, scale, 2); 5543c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale); 5553c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale); 55662c23b28SMark if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 55762c23b28SMark PetscFunctionReturn(0); 55862c23b28SMark } 55962c23b28SMark 56062c23b28SMark /*@ 56162c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 56262c23b28SMark 56362c23b28SMark Not Collective 56462c23b28SMark 5654165533cSJose E. Roman Input Parameter: 56662c23b28SMark . adapt - adaptive controller context 56762c23b28SMark 5684165533cSJose E. Roman Output Parameter: 569ee300463SSatish Balay . scale - scale factor 57062c23b28SMark 57162c23b28SMark Level: intermediate 57262c23b28SMark 573db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 57462c23b28SMark @*/ 5759371c9d4SSatish Balay PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale) { 57662c23b28SMark PetscFunctionBegin; 57762c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 57862c23b28SMark if (scale) PetscValidRealPointer(scale, 2); 57962c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 58062c23b28SMark PetscFunctionReturn(0); 58162c23b28SMark } 58262c23b28SMark 58362c23b28SMark /*@ 5841917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 5851917a363SLisandro Dalcin 5861917a363SLisandro Dalcin Logically collective on TSAdapt 5871c3436cfSJed Brown 5884165533cSJose E. Roman Input Parameters: 589552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5901c3436cfSJed Brown . hmin - minimum time step 5911c3436cfSJed Brown - hmax - maximum time step 5921c3436cfSJed Brown 5931c3436cfSJed Brown Options Database Keys: 594ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 595ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 5961c3436cfSJed Brown 5971c3436cfSJed Brown Level: intermediate 5981c3436cfSJed Brown 599db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 6001c3436cfSJed Brown @*/ 6019371c9d4SSatish Balay PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax) { 6021c3436cfSJed Brown PetscFunctionBegin; 6034782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 604b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmin, 2); 605b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmax, 3); 6063c633725SBarry Smith PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin); 6073c633725SBarry Smith PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax); 608b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 609b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 610b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 611b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6123c633725SBarry 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); 6131917a363SLisandro Dalcin PetscFunctionReturn(0); 6141917a363SLisandro Dalcin } 6151917a363SLisandro Dalcin 6161917a363SLisandro Dalcin /*@ 6171917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 6181917a363SLisandro Dalcin 6191917a363SLisandro Dalcin Not Collective 6201917a363SLisandro Dalcin 6214165533cSJose E. Roman Input Parameter: 6221917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6231917a363SLisandro Dalcin 6244165533cSJose E. Roman Output Parameters: 6251917a363SLisandro Dalcin + hmin - minimum time step 6261917a363SLisandro Dalcin - hmax - maximum time step 6271917a363SLisandro Dalcin 6281917a363SLisandro Dalcin Level: intermediate 6291917a363SLisandro Dalcin 630db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 6311917a363SLisandro Dalcin @*/ 6329371c9d4SSatish Balay PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax) { 6331917a363SLisandro Dalcin PetscFunctionBegin; 6341917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 6351917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin, 2); 6361917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax, 3); 6371917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6381917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6391c3436cfSJed Brown PetscFunctionReturn(0); 6401c3436cfSJed Brown } 6411c3436cfSJed Brown 642e55864a3SBarry Smith /* 64384df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 64484df9cb4SJed Brown 64584df9cb4SJed Brown Collective on TSAdapt 64684df9cb4SJed Brown 64784df9cb4SJed Brown Input Parameter: 64884df9cb4SJed Brown . adapt - the TSAdapt context 64984df9cb4SJed Brown 65084df9cb4SJed Brown Options Database Keys: 6511917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 652bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6531917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6541917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6551917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 65623de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 65723de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 65823de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 659de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 660de50f1caSBarry 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 66184df9cb4SJed Brown 66284df9cb4SJed Brown Level: advanced 66384df9cb4SJed Brown 66484df9cb4SJed Brown Notes: 66584df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 66684df9cb4SJed Brown 667db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 668db781477SPatrick Sanan `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 669e55864a3SBarry Smith */ 6709371c9d4SSatish Balay PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject) { 67184df9cb4SJed Brown char type[256] = TSADAPTBASIC; 67262c23b28SMark PetscReal safety, reject_safety, clip[2], scale, hmin, hmax; 6731c3436cfSJed Brown PetscBool set, flg; 6741917a363SLisandro Dalcin PetscInt two; 67584df9cb4SJed Brown 67684df9cb4SJed Brown PetscFunctionBegin; 677dbbe0bcdSBarry Smith PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 67884df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6791566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 680d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options"); 6819566063dSJacob 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)); 682*48a46eb9SPierre Jolivet if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type)); 6831917a363SLisandro Dalcin 6849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set)); 6859566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg)); 6861917a363SLisandro Dalcin 6879371c9d4SSatish Balay safety = adapt->safety; 6889371c9d4SSatish Balay reject_safety = adapt->reject_safety; 6899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set)); 6909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg)); 6919566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety)); 6921917a363SLisandro Dalcin 6939371c9d4SSatish Balay two = 2; 6949371c9d4SSatish Balay clip[0] = adapt->clip[0]; 6959371c9d4SSatish Balay clip[1] = adapt->clip[1]; 6969566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set)); 6973c633725SBarry Smith PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip"); 6989566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1])); 6991917a363SLisandro Dalcin 7009371c9d4SSatish Balay hmin = adapt->dt_min; 7019371c9d4SSatish Balay hmax = adapt->dt_max; 7029566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set)); 7039566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg)); 7049566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax)); 7051917a363SLisandro Dalcin 7069566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set)); 7079566063dSJacob 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)); 708d580f011SEmil Constantinescu 7099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set)); 7109566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale)); 7111917a363SLisandro Dalcin 7129566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL)); 7133c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported"); 7141917a363SLisandro Dalcin 7159566063dSJacob 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)); 716de50f1caSBarry Smith 7179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set)); 7189566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt, flg)); 7191917a363SLisandro Dalcin 720dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 721d0609cedSBarry Smith PetscOptionsHeadEnd(); 72284df9cb4SJed Brown PetscFunctionReturn(0); 72384df9cb4SJed Brown } 72484df9cb4SJed Brown 72584df9cb4SJed Brown /*@ 72684df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 72784df9cb4SJed Brown 7281917a363SLisandro Dalcin Logically collective on TSAdapt 72984df9cb4SJed Brown 7304165533cSJose E. Roman Input Parameter: 73184df9cb4SJed Brown . adapt - adaptive controller 73284df9cb4SJed Brown 73384df9cb4SJed Brown Level: developer 73484df9cb4SJed Brown 735db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 73684df9cb4SJed Brown @*/ 7379371c9d4SSatish Balay PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) { 73884df9cb4SJed Brown PetscFunctionBegin; 7394782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 7409566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates))); 74184df9cb4SJed Brown PetscFunctionReturn(0); 74284df9cb4SJed Brown } 74384df9cb4SJed Brown 74484df9cb4SJed Brown /*@C 74584df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 74684df9cb4SJed Brown 7471917a363SLisandro Dalcin Logically collective on TSAdapt 74884df9cb4SJed Brown 7494165533cSJose E. Roman Input Parameters: 750552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 75184df9cb4SJed Brown . name - name of the candidate scheme to add 75284df9cb4SJed Brown . order - order of the candidate scheme 75384df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7548d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 75584df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 75684df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 75784df9cb4SJed Brown 75884df9cb4SJed Brown Note: 75984df9cb4SJed Brown This routine is not available in Fortran. 76084df9cb4SJed Brown 76184df9cb4SJed Brown Level: developer 76284df9cb4SJed Brown 763db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 76484df9cb4SJed Brown @*/ 7659371c9d4SSatish Balay PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse) { 76684df9cb4SJed Brown PetscInt c; 76784df9cb4SJed Brown 76884df9cb4SJed Brown PetscFunctionBegin; 76984df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 77063a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order); 77184df9cb4SJed Brown if (inuse) { 7723c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 77384df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 77484df9cb4SJed Brown } 7751c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7761c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 777bbd56ea5SKarl Rupp 77884df9cb4SJed Brown adapt->candidates.name[c] = name; 77984df9cb4SJed Brown adapt->candidates.order[c] = order; 78084df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7818d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 78284df9cb4SJed Brown adapt->candidates.cost[c] = cost; 78384df9cb4SJed Brown adapt->candidates.n++; 78484df9cb4SJed Brown PetscFunctionReturn(0); 78584df9cb4SJed Brown } 78684df9cb4SJed Brown 7878d59e960SJed Brown /*@C 7888d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 7898d59e960SJed Brown 7908d59e960SJed Brown Not Collective 7918d59e960SJed Brown 7924165533cSJose E. Roman Input Parameter: 7938d59e960SJed Brown . adapt - time step adaptivity context 7948d59e960SJed Brown 7954165533cSJose E. Roman Output Parameters: 7968d59e960SJed Brown + n - number of candidate schemes, always at least 1 7978d59e960SJed Brown . order - the order of each candidate scheme 7988d59e960SJed Brown . stageorder - the stage order of each candidate scheme 7998d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8008d59e960SJed Brown - cost - the relative cost of each scheme 8018d59e960SJed Brown 8028d59e960SJed Brown Level: developer 8038d59e960SJed Brown 8048d59e960SJed Brown Note: 8058d59e960SJed Brown The current scheme is always returned in the first slot 8068d59e960SJed Brown 807db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 8088d59e960SJed Brown @*/ 8099371c9d4SSatish Balay PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost) { 8108d59e960SJed Brown PetscFunctionBegin; 8118d59e960SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 8128d59e960SJed Brown if (n) *n = adapt->candidates.n; 8138d59e960SJed Brown if (order) *order = adapt->candidates.order; 8148d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8158d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8168d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8178d59e960SJed Brown PetscFunctionReturn(0); 8188d59e960SJed Brown } 8198d59e960SJed Brown 82084df9cb4SJed Brown /*@C 82184df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 82284df9cb4SJed Brown 8231917a363SLisandro Dalcin Collective on TSAdapt 82484df9cb4SJed Brown 8254165533cSJose E. Roman Input Parameters: 82684df9cb4SJed Brown + adapt - adaptive contoller 82797bb3fdcSJose E. Roman . ts - time stepper 82884df9cb4SJed Brown - h - current step size 82984df9cb4SJed Brown 8304165533cSJose E. Roman Output Parameters: 8311566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 83284df9cb4SJed Brown . next_h - step size to use for the next step 83384df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 83484df9cb4SJed Brown 8351c3436cfSJed Brown Note: 8361c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 8371c3436cfSJed Brown being retried after an initial rejection. 8381c3436cfSJed Brown 83984df9cb4SJed Brown Level: developer 84084df9cb4SJed Brown 841db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 84284df9cb4SJed Brown @*/ 8439371c9d4SSatish Balay PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept) { 8441566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8451566a47fSLisandro Dalcin PetscInt scheme = 0; 8460b99f514SJed Brown PetscReal wlte = -1.0; 8475e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8485e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 84984df9cb4SJed Brown 85084df9cb4SJed Brown PetscFunctionBegin; 85184df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 85284df9cb4SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 8531566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc, 4); 854dadcf809SJacob Faibussowitsch PetscValidRealPointer(next_h, 5); 855064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept, 6); 8561566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8571566a47fSLisandro Dalcin 8581566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8591566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8601566a47fSLisandro Dalcin *next_h = h; 8611566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8621566a47fSLisandro Dalcin PetscFunctionReturn(0); 8631566a47fSLisandro Dalcin } 8641566a47fSLisandro Dalcin 865dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter); 86663a3b9bcSJacob 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); 8673c633725SBarry Smith PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h); 8681566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8691566a47fSLisandro Dalcin 87049354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 87136b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 87236b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 8734a658b32SHong Zhang PetscReal tend = t + h, tmax, hmax; 874fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 875fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 8764a658b32SHong Zhang 8774a658b32SHong Zhang if (ts->tspan) { 878e1db57b0SHong Zhang if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */ 8794a658b32SHong Zhang if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1]; 8804a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */ 8814a658b32SHong Zhang else tmax = ts->tspan->span_times[ts->tspan->spanctr]; 8824a658b32SHong Zhang } else tmax = ts->max_time; 8834a658b32SHong Zhang hmax = tmax - t; 88436b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 885fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2; 88636b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax; 887e1db57b0SHong Zhang /* if step size is changed to match a span time point */ 888e1db57b0SHong Zhang if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h; 889e1db57b0SHong Zhang /* reset time step after a span time point */ 890e1db57b0SHong Zhang if (ts->tspan && h == *next_h && adapt->dt_span_cached && PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) { 891e1db57b0SHong Zhang *next_h = adapt->dt_span_cached; 892e1db57b0SHong Zhang adapt->dt_span_cached = 0; 89349354f04SShri Abhyankar } 894e1db57b0SHong Zhang } 8951c3436cfSJed Brown if (adapt->monitor) { 8961566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 8980b99f514SJed Brown if (wlte < 0) { 8999371c9d4SSatish 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", 9009371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h)); 9010b99f514SJed Brown } else { 9029371c9d4SSatish 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", 9039371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter)); 9040b99f514SJed Brown } 9059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9061c3436cfSJed Brown } 90784df9cb4SJed Brown PetscFunctionReturn(0); 90884df9cb4SJed Brown } 90984df9cb4SJed Brown 91097335746SJed Brown /*@ 911de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 912de50f1caSBarry Smith before increasing the time step. 913de50f1caSBarry Smith 914de50f1caSBarry Smith Logicially Collective on TSAdapt 915de50f1caSBarry Smith 9164165533cSJose E. Roman Input Parameters: 917de50f1caSBarry Smith + adapt - adaptive controller context 918de50f1caSBarry Smith - cnt - the number of timesteps 919de50f1caSBarry Smith 920de50f1caSBarry Smith Options Database Key: 921de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 922de50f1caSBarry Smith 923de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 924de50f1caSBarry Smith The successful use of this option is problem dependent 925de50f1caSBarry Smith 926de50f1caSBarry Smith Developer Note: there is no theory to support this option 927de50f1caSBarry Smith 928de50f1caSBarry Smith Level: advanced 929de50f1caSBarry Smith 930de50f1caSBarry Smith .seealso: 931de50f1caSBarry Smith @*/ 9329371c9d4SSatish Balay PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt) { 933de50f1caSBarry Smith PetscFunctionBegin; 934de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 935de50f1caSBarry Smith PetscFunctionReturn(0); 936de50f1caSBarry Smith } 937de50f1caSBarry Smith 938de50f1caSBarry Smith /*@ 9396bc98fa9SBarry 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) 94097335746SJed Brown 9411917a363SLisandro Dalcin Collective on TSAdapt 94297335746SJed Brown 9434165533cSJose E. Roman Input Parameters: 94497335746SJed Brown + adapt - adaptive controller context 945b295832fSPierre Barbier de Reuille . ts - time stepper 946b295832fSPierre Barbier de Reuille . t - Current simulation time 947b295832fSPierre Barbier de Reuille - Y - Current solution vector 94897335746SJed Brown 9494165533cSJose E. Roman Output Parameter: 95097335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 95197335746SJed Brown 95297335746SJed Brown Level: developer 95397335746SJed Brown 95497335746SJed Brown .seealso: 95597335746SJed Brown @*/ 9569371c9d4SSatish Balay PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) { 9571566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 95897335746SJed Brown 95997335746SJed Brown PetscFunctionBegin; 9604782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 9614782b174SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 962064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept, 5); 9631566a47fSLisandro Dalcin 9649566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 96597335746SJed Brown if (snesreason < 0) { 96697335746SJed Brown *accept = PETSC_FALSE; 9676de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 96897335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 96963a3b9bcSJacob 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)); 97097335746SJed Brown if (adapt->monitor) { 9719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9729371c9d4SSatish 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, 9739371c9d4SSatish Balay (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures)); 9749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 97597335746SJed Brown } 976cb9d8021SPierre Barbier de Reuille } 977cb9d8021SPierre Barbier de Reuille } else { 9781566a47fSLisandro Dalcin *accept = PETSC_TRUE; 9799566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts, t, Y, accept)); 980cb9d8021SPierre Barbier de Reuille if (*accept && adapt->checkstage) { 9819566063dSJacob Faibussowitsch PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept)); 9826bc98fa9SBarry Smith if (!*accept) { 98363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n", ts->steps)); 9846bc98fa9SBarry Smith if (adapt->monitor) { 9859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 98663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected by user function provided by TSSetFunctionDomainError()\n", ((PetscObject)adapt)->type_name, ts->steps)); 9879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 9886bc98fa9SBarry Smith } 9896bc98fa9SBarry Smith } 990cb9d8021SPierre Barbier de Reuille } 991cb9d8021SPierre Barbier de Reuille } 992cb9d8021SPierre Barbier de Reuille 9931566a47fSLisandro Dalcin if (!(*accept) && !ts->reason) { 9941566a47fSLisandro Dalcin PetscReal dt, new_dt; 9959566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 996cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 9979566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, new_dt)); 998de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 99997335746SJed Brown if (adapt->monitor) { 10009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 100163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ts->steps, SNESConvergedReasons[snesreason], (double)ts->ptime, (double)dt, (double)new_dt)); 10029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 100397335746SJed Brown } 100497335746SJed Brown } 100597335746SJed Brown PetscFunctionReturn(0); 100697335746SJed Brown } 100797335746SJed Brown 100884df9cb4SJed Brown /*@ 100984df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 101084df9cb4SJed Brown 1011d083f849SBarry Smith Collective 101284df9cb4SJed Brown 101384df9cb4SJed Brown Input Parameter: 101484df9cb4SJed Brown . comm - The communicator 101584df9cb4SJed Brown 101684df9cb4SJed Brown Output Parameter: 101784df9cb4SJed Brown . adapt - new TSAdapt object 101884df9cb4SJed Brown 101984df9cb4SJed Brown Level: developer 102084df9cb4SJed Brown 102184df9cb4SJed Brown Notes: 102284df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 102384df9cb4SJed Brown 1024db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 102584df9cb4SJed Brown @*/ 10269371c9d4SSatish Balay PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt) { 102784df9cb4SJed Brown TSAdapt adapt; 102884df9cb4SJed Brown 102984df9cb4SJed Brown PetscFunctionBegin; 1030064a246eSJacob Faibussowitsch PetscValidPointer(inadapt, 2); 10313b3bcf4cSLisandro Dalcin *inadapt = NULL; 10329566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 10333b3bcf4cSLisandro Dalcin 10349566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView)); 10351c3436cfSJed Brown 1036bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 10371917a363SLisandro Dalcin adapt->safety = 0.9; 10381917a363SLisandro Dalcin adapt->reject_safety = 0.5; 10391917a363SLisandro Dalcin adapt->clip[0] = 0.1; 10401917a363SLisandro Dalcin adapt->clip[1] = 10.; 10411c3436cfSJed Brown adapt->dt_min = 1e-20; 10421917a363SLisandro Dalcin adapt->dt_max = 1e+20; 10431c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1044d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 104597335746SJed Brown adapt->scale_solve_failed = 0.25; 1046fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1047fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1048fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1049fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 10507619abb3SShri adapt->wnormtype = NORM_2; 1051de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 10521917a363SLisandro Dalcin 105384df9cb4SJed Brown *inadapt = adapt; 105484df9cb4SJed Brown PetscFunctionReturn(0); 105584df9cb4SJed Brown } 1056