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 @*/ 43bdf89e91SBarry Smith PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 4484df9cb4SJed Brown { 4584df9cb4SJed Brown PetscFunctionBegin; 469566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 479566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSAdaptList,sname,function)); 4884df9cb4SJed Brown PetscFunctionReturn(0); 4984df9cb4SJed Brown } 5084df9cb4SJed Brown 5184df9cb4SJed Brown /*@C 5284df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5384df9cb4SJed Brown 5484df9cb4SJed Brown Not Collective 5584df9cb4SJed Brown 5684df9cb4SJed Brown Level: advanced 5784df9cb4SJed Brown 58db781477SPatrick Sanan .seealso: `TSAdaptRegisterDestroy()` 5984df9cb4SJed Brown @*/ 60607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6184df9cb4SJed Brown { 6284df9cb4SJed Brown PetscFunctionBegin; 630f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 640f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 659566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None)); 669566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic)); 679566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP)); 689566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL)); 699566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE)); 709566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTHISTORY,TSAdaptCreate_History)); 7184df9cb4SJed Brown PetscFunctionReturn(0); 7284df9cb4SJed Brown } 7384df9cb4SJed Brown 7484df9cb4SJed Brown /*@C 75560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 7684df9cb4SJed Brown called from PetscFinalize(). 7784df9cb4SJed Brown 7884df9cb4SJed Brown Level: developer 7984df9cb4SJed Brown 80db781477SPatrick Sanan .seealso: `PetscFinalize()` 8184df9cb4SJed Brown @*/ 8284df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8384df9cb4SJed Brown { 8484df9cb4SJed Brown PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSAdaptList)); 8684df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 8784df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 8884df9cb4SJed Brown PetscFunctionReturn(0); 8984df9cb4SJed Brown } 9084df9cb4SJed Brown 9184df9cb4SJed Brown /*@C 9284df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 938a690491SBarry Smith called from TSInitializePackage(). 9484df9cb4SJed Brown 9584df9cb4SJed Brown Level: developer 9684df9cb4SJed Brown 97db781477SPatrick Sanan .seealso: `PetscInitialize()` 9884df9cb4SJed Brown @*/ 99607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 10084df9cb4SJed Brown { 10184df9cb4SJed Brown PetscFunctionBegin; 10284df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 10384df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 1049566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID)); 1059566063dSJacob Faibussowitsch PetscCall(TSAdaptRegisterAll()); 1069566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage)); 10784df9cb4SJed Brown PetscFunctionReturn(0); 10884df9cb4SJed Brown } 10984df9cb4SJed Brown 1107eef6055SBarry Smith /*@C 1117eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1127eef6055SBarry Smith 1137eef6055SBarry Smith Logicially Collective on TSAdapt 1147eef6055SBarry Smith 115d8d19677SJose E. Roman Input Parameters: 116d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1177eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1187eef6055SBarry Smith 1197eef6055SBarry Smith Options Database: 120ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1217eef6055SBarry Smith 1227eef6055SBarry Smith Level: intermediate 1237eef6055SBarry Smith 124db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()` 1257eef6055SBarry Smith @*/ 12619fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 12784df9cb4SJed Brown { 128ef749922SLisandro Dalcin PetscBool match; 1295f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSAdapt); 13084df9cb4SJed Brown 13184df9cb4SJed Brown PetscFunctionBegin; 1324782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 133b92453a8SLisandro Dalcin PetscValidCharPointer(type,2); 1349566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt,type,&match)); 135ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1369566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSAdaptList,type,&r)); 1373c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 1389566063dSJacob Faibussowitsch if (adapt->ops->destroy) PetscCall((*adapt->ops->destroy)(adapt)); 1399566063dSJacob Faibussowitsch PetscCall(PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps))); 1409566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt,type)); 1419566063dSJacob Faibussowitsch PetscCall((*r)(adapt)); 14284df9cb4SJed Brown PetscFunctionReturn(0); 14384df9cb4SJed Brown } 14484df9cb4SJed Brown 145d0288e4fSLisandro Dalcin /*@C 146d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 147d0288e4fSLisandro Dalcin 148d0288e4fSLisandro Dalcin Not Collective 149d0288e4fSLisandro Dalcin 150d0288e4fSLisandro Dalcin Input Parameter: 151d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 152d0288e4fSLisandro Dalcin 153d0288e4fSLisandro Dalcin Output Parameter: 154d0288e4fSLisandro Dalcin . type - The name of TS adapter method 155d0288e4fSLisandro Dalcin 156d0288e4fSLisandro Dalcin Level: intermediate 157d0288e4fSLisandro Dalcin 158db781477SPatrick Sanan .seealso `TSAdaptSetType()` 159d0288e4fSLisandro Dalcin @*/ 160d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 161d0288e4fSLisandro Dalcin { 162d0288e4fSLisandro Dalcin PetscFunctionBegin; 163d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 164d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 165d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 166d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 167d0288e4fSLisandro Dalcin } 168d0288e4fSLisandro Dalcin 16984df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 17084df9cb4SJed Brown { 17184df9cb4SJed Brown PetscFunctionBegin; 1724782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix)); 17484df9cb4SJed Brown PetscFunctionReturn(0); 17584df9cb4SJed Brown } 17684df9cb4SJed Brown 177ad6bc421SBarry Smith /*@C 178ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 179ad6bc421SBarry Smith 180ad6bc421SBarry Smith Collective on PetscViewer 181ad6bc421SBarry Smith 182ad6bc421SBarry Smith Input Parameters: 183ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 184ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 185ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 186ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 187ad6bc421SBarry Smith 188ad6bc421SBarry Smith Level: intermediate 189ad6bc421SBarry Smith 190ad6bc421SBarry Smith Notes: 191ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 192ad6bc421SBarry Smith 193ad6bc421SBarry Smith Notes for advanced users: 194ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 195ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 196ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 197ad6bc421SBarry Smith format is 198ad6bc421SBarry Smith .vb 199ad6bc421SBarry Smith has not yet been determined 200ad6bc421SBarry Smith .ve 201ad6bc421SBarry Smith 202db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()` 203ad6bc421SBarry Smith @*/ 2044782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 205ad6bc421SBarry Smith { 206ad6bc421SBarry Smith PetscBool isbinary; 207ad6bc421SBarry Smith char type[256]; 208ad6bc421SBarry Smith 209ad6bc421SBarry Smith PetscFunctionBegin; 2104782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 211ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 2133c633725SBarry Smith PetscCheck(isbinary,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 214ad6bc421SBarry Smith 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR)); 2169566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt,type)); 2179566063dSJacob Faibussowitsch if (adapt->ops->load) PetscCall((*adapt->ops->load)(adapt,viewer)); 218ad6bc421SBarry Smith PetscFunctionReturn(0); 219ad6bc421SBarry Smith } 220ad6bc421SBarry Smith 22184df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 22284df9cb4SJed Brown { 2231c167fc2SEmil Constantinescu PetscBool iascii,isbinary,isnone,isglee; 22484df9cb4SJed Brown 22584df9cb4SJed Brown PetscFunctionBegin; 2264782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2279566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer)); 2284782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2294782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 23284df9cb4SJed Brown if (iascii) { 2339566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer)); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt,TSADAPTGLEE,&isglee)); 2361917a363SLisandro Dalcin if (!isnone) { 2379566063dSJacob Faibussowitsch if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer," always accepting steps\n")); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety)); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety)); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1])); 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0])); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min)); 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," maximum solution absolute value to be ignored %g\n",(double)adapt->ignore_max)); 2451c167fc2SEmil Constantinescu } 2461c167fc2SEmil Constantinescu if (isglee) { 2471c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," GLEE uses local error control\n")); 2491c167fc2SEmil Constantinescu } else { 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," GLEE uses global error control\n")); 2511c167fc2SEmil Constantinescu } 2521917a363SLisandro Dalcin } 25384df9cb4SJed Brown if (adapt->ops->view) { 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2559566063dSJacob Faibussowitsch PetscCall((*adapt->ops->view)(adapt,viewer)); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 25784df9cb4SJed Brown } 258ad6bc421SBarry Smith } else if (isbinary) { 259ad6bc421SBarry Smith char type[256]; 260ad6bc421SBarry Smith 261ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 2629566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type,((PetscObject)adapt)->type_name,256)); 2639566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR)); 2641baa6e33SBarry Smith } else if (adapt->ops->view) PetscCall((*adapt->ops->view)(adapt,viewer)); 26584df9cb4SJed Brown PetscFunctionReturn(0); 26684df9cb4SJed Brown } 26784df9cb4SJed Brown 268099af0a3SLisandro Dalcin /*@ 269099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 270099af0a3SLisandro Dalcin 271099af0a3SLisandro Dalcin Collective on TS 272099af0a3SLisandro Dalcin 273099af0a3SLisandro Dalcin Input Parameter: 274099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 275099af0a3SLisandro Dalcin 276099af0a3SLisandro Dalcin Level: developer 277099af0a3SLisandro Dalcin 278db781477SPatrick Sanan .seealso: `TSAdaptCreate()`, `TSAdaptDestroy()` 279099af0a3SLisandro Dalcin @*/ 280099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 281099af0a3SLisandro Dalcin { 282099af0a3SLisandro Dalcin PetscFunctionBegin; 283099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2849566063dSJacob Faibussowitsch if (adapt->ops->reset) PetscCall((*adapt->ops->reset)(adapt)); 285099af0a3SLisandro Dalcin PetscFunctionReturn(0); 286099af0a3SLisandro Dalcin } 287099af0a3SLisandro Dalcin 28884df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 28984df9cb4SJed Brown { 29084df9cb4SJed Brown PetscFunctionBegin; 29184df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 29284df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 2934782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 294099af0a3SLisandro Dalcin 2959566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt)); 296099af0a3SLisandro Dalcin 2979566063dSJacob Faibussowitsch if ((*adapt)->ops->destroy) PetscCall((*(*adapt)->ops->destroy)(*adapt)); 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 2999566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 30084df9cb4SJed Brown PetscFunctionReturn(0); 30184df9cb4SJed Brown } 30284df9cb4SJed Brown 3031c3436cfSJed Brown /*@ 3041c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3051c3436cfSJed Brown 3061c3436cfSJed Brown Collective on TSAdapt 3071c3436cfSJed Brown 3084165533cSJose E. Roman Input Parameters: 3091c3436cfSJed Brown + adapt - adaptive controller context 3101c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3111c3436cfSJed Brown 312bf997491SLisandro Dalcin Options Database Keys: 313ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 314bf997491SLisandro Dalcin 3151c3436cfSJed Brown Level: intermediate 3161c3436cfSJed Brown 317db781477SPatrick Sanan .seealso: `TSAdaptChoose()` 3181c3436cfSJed Brown @*/ 3191c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3201c3436cfSJed Brown { 3211c3436cfSJed Brown PetscFunctionBegin; 3224782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3234782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3241c3436cfSJed Brown if (flg) { 3259566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor)); 3261c3436cfSJed Brown } else { 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor)); 3281c3436cfSJed Brown } 3291c3436cfSJed Brown PetscFunctionReturn(0); 3301c3436cfSJed Brown } 3311c3436cfSJed Brown 3320873808bSJed Brown /*@C 333bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3340873808bSJed Brown 3350873808bSJed Brown Logically collective on TSAdapt 3360873808bSJed Brown 3374165533cSJose E. Roman Input Parameters: 3380873808bSJed Brown + adapt - adaptive controller context 3390873808bSJed Brown - func - stage check function 3400873808bSJed Brown 3410873808bSJed Brown Arguments of func: 3420873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3430873808bSJed Brown 3440873808bSJed Brown + adapt - adaptive controller context 3450873808bSJed Brown . ts - time stepping context 3460873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3470873808bSJed Brown 3480873808bSJed Brown Level: advanced 3490873808bSJed Brown 350db781477SPatrick Sanan .seealso: `TSAdaptChoose()` 3510873808bSJed Brown @*/ 352b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3530873808bSJed Brown { 3540873808bSJed Brown PetscFunctionBegin; 3550873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 35668ae941aSLisandro Dalcin adapt->checkstage = func; 3570873808bSJed Brown PetscFunctionReturn(0); 3580873808bSJed Brown } 3590873808bSJed Brown 3601c3436cfSJed Brown /*@ 361bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 362bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 363bf997491SLisandro Dalcin 3641917a363SLisandro Dalcin Logically collective on TSAdapt 365bf997491SLisandro Dalcin 3664165533cSJose E. Roman Input Parameters: 367bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 368bf997491SLisandro Dalcin - flag - whether to always accept steps 369bf997491SLisandro Dalcin 370bf997491SLisandro Dalcin Options Database Keys: 371ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 372bf997491SLisandro Dalcin 373bf997491SLisandro Dalcin Level: intermediate 374bf997491SLisandro Dalcin 375db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptChoose()` 376bf997491SLisandro Dalcin @*/ 377bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 378bf997491SLisandro Dalcin { 379bf997491SLisandro Dalcin PetscFunctionBegin; 380bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 381bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 382bf997491SLisandro Dalcin adapt->always_accept = flag; 383bf997491SLisandro Dalcin PetscFunctionReturn(0); 384bf997491SLisandro Dalcin } 385bf997491SLisandro Dalcin 386bf997491SLisandro Dalcin /*@ 3871917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 3881c3436cfSJed Brown 3891917a363SLisandro Dalcin Logically collective on TSAdapt 3901917a363SLisandro Dalcin 3914165533cSJose E. Roman Input Parameters: 3921917a363SLisandro Dalcin + adapt - adaptive controller context 3931917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 394ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 3951917a363SLisandro Dalcin 3961917a363SLisandro Dalcin Options Database Keys: 397ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 398ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 3991917a363SLisandro Dalcin 4001917a363SLisandro Dalcin Level: intermediate 4011917a363SLisandro Dalcin 402db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 4031917a363SLisandro Dalcin @*/ 4041917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4051917a363SLisandro Dalcin { 4061917a363SLisandro Dalcin PetscFunctionBegin; 4071917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4081917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4091917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4103c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4113c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety); 4123c633725SBarry 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); 4133c633725SBarry 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); 4141917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4151917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4161917a363SLisandro Dalcin PetscFunctionReturn(0); 4171917a363SLisandro Dalcin } 4181917a363SLisandro Dalcin 4191917a363SLisandro Dalcin /*@ 4201917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4211917a363SLisandro Dalcin 4221917a363SLisandro Dalcin Not Collective 4231917a363SLisandro Dalcin 4244165533cSJose E. Roman Input Parameter: 4251917a363SLisandro Dalcin . adapt - adaptive controller context 4261917a363SLisandro Dalcin 4274165533cSJose E. Roman Output Parameters: 4281917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4291917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4301917a363SLisandro Dalcin 4311917a363SLisandro Dalcin Level: intermediate 4321917a363SLisandro Dalcin 433db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 4341917a363SLisandro Dalcin @*/ 4351917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4361917a363SLisandro Dalcin { 4371917a363SLisandro Dalcin PetscFunctionBegin; 4381917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4391917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4401917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4411917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4421917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4431917a363SLisandro Dalcin PetscFunctionReturn(0); 4441917a363SLisandro Dalcin } 4451917a363SLisandro Dalcin 4461917a363SLisandro Dalcin /*@ 44776cddca1SEmil 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. 44876cddca1SEmil Constantinescu 44976cddca1SEmil Constantinescu Logically collective on TSAdapt 45076cddca1SEmil Constantinescu 4514165533cSJose E. Roman Input Parameters: 45276cddca1SEmil Constantinescu + adapt - adaptive controller context 45376cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 45476cddca1SEmil Constantinescu 45576cddca1SEmil Constantinescu Options Database Keys: 45676cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 45776cddca1SEmil Constantinescu 45876cddca1SEmil Constantinescu Level: intermediate 45976cddca1SEmil Constantinescu 460db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 46176cddca1SEmil Constantinescu @*/ 46276cddca1SEmil Constantinescu PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore) 46376cddca1SEmil Constantinescu { 46476cddca1SEmil Constantinescu PetscFunctionBegin; 46576cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 46676cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt,max_ignore,2); 46776cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 46876cddca1SEmil Constantinescu PetscFunctionReturn(0); 46976cddca1SEmil Constantinescu } 47076cddca1SEmil Constantinescu 47176cddca1SEmil Constantinescu /*@ 47276cddca1SEmil 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). 47376cddca1SEmil Constantinescu 47476cddca1SEmil Constantinescu Not Collective 47576cddca1SEmil Constantinescu 4764165533cSJose E. Roman Input Parameter: 47776cddca1SEmil Constantinescu . adapt - adaptive controller context 47876cddca1SEmil Constantinescu 4794165533cSJose E. Roman Output Parameter: 48076cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 48176cddca1SEmil Constantinescu 48276cddca1SEmil Constantinescu Level: intermediate 48376cddca1SEmil Constantinescu 484db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 48576cddca1SEmil Constantinescu @*/ 48676cddca1SEmil Constantinescu PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore) 48776cddca1SEmil Constantinescu { 48876cddca1SEmil Constantinescu PetscFunctionBegin; 48976cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 49076cddca1SEmil Constantinescu PetscValidRealPointer(max_ignore,2); 49176cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 49276cddca1SEmil Constantinescu PetscFunctionReturn(0); 49376cddca1SEmil Constantinescu } 49476cddca1SEmil Constantinescu 49576cddca1SEmil Constantinescu /*@ 4961917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4971917a363SLisandro Dalcin 4981917a363SLisandro Dalcin Logically collective on TSAdapt 4991917a363SLisandro Dalcin 5004165533cSJose E. Roman Input Parameters: 5011917a363SLisandro Dalcin + adapt - adaptive controller context 5021917a363SLisandro Dalcin . low - admissible decrease factor 503ec18b7bcSLisandro Dalcin - high - admissible increase factor 5041917a363SLisandro Dalcin 5051917a363SLisandro Dalcin Options Database Keys: 506ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5071917a363SLisandro Dalcin 5081917a363SLisandro Dalcin Level: intermediate 5091917a363SLisandro Dalcin 510db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 5111917a363SLisandro Dalcin @*/ 5121917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 5131917a363SLisandro Dalcin { 5141917a363SLisandro Dalcin PetscFunctionBegin; 5151917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5161917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 5171917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 5183c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 5193c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 5203c633725SBarry Smith PetscCheck(high == PETSC_DEFAULT || high >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be greater than one",(double)high); 5211917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 5221917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 5231917a363SLisandro Dalcin PetscFunctionReturn(0); 5241917a363SLisandro Dalcin } 5251917a363SLisandro Dalcin 5261917a363SLisandro Dalcin /*@ 5271917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5281917a363SLisandro Dalcin 5291917a363SLisandro Dalcin Not Collective 5301917a363SLisandro Dalcin 5314165533cSJose E. Roman Input Parameter: 5321917a363SLisandro Dalcin . adapt - adaptive controller context 5331917a363SLisandro Dalcin 5344165533cSJose E. Roman Output Parameters: 5351917a363SLisandro Dalcin + low - optional, admissible decrease factor 5361917a363SLisandro Dalcin - high - optional, admissible increase factor 5371917a363SLisandro Dalcin 5381917a363SLisandro Dalcin Level: intermediate 5391917a363SLisandro Dalcin 540db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 5411917a363SLisandro Dalcin @*/ 5421917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5431917a363SLisandro Dalcin { 5441917a363SLisandro Dalcin PetscFunctionBegin; 5451917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5461917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5471917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5481917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5491917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5501917a363SLisandro Dalcin PetscFunctionReturn(0); 5511917a363SLisandro Dalcin } 5521917a363SLisandro Dalcin 5531917a363SLisandro Dalcin /*@ 55462c23b28SMark TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 55562c23b28SMark 55662c23b28SMark Logically collective on TSAdapt 55762c23b28SMark 5584165533cSJose E. Roman Input Parameters: 55962c23b28SMark + adapt - adaptive controller context 560ee300463SSatish Balay - scale - scale 56162c23b28SMark 56262c23b28SMark Options Database Keys: 56362c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 56462c23b28SMark 56562c23b28SMark Level: intermediate 56662c23b28SMark 567db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 56862c23b28SMark @*/ 56962c23b28SMark PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt,PetscReal scale) 57062c23b28SMark { 57162c23b28SMark PetscFunctionBegin; 57262c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 57362c23b28SMark PetscValidLogicalCollectiveReal(adapt,scale,2); 5743c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be positive",(double)scale); 5753c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be less than one",(double)scale); 57662c23b28SMark if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 57762c23b28SMark PetscFunctionReturn(0); 57862c23b28SMark } 57962c23b28SMark 58062c23b28SMark /*@ 58162c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 58262c23b28SMark 58362c23b28SMark Not Collective 58462c23b28SMark 5854165533cSJose E. Roman Input Parameter: 58662c23b28SMark . adapt - adaptive controller context 58762c23b28SMark 5884165533cSJose E. Roman Output Parameter: 589ee300463SSatish Balay . scale - scale factor 59062c23b28SMark 59162c23b28SMark Level: intermediate 59262c23b28SMark 593db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 59462c23b28SMark @*/ 59562c23b28SMark PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt,PetscReal *scale) 59662c23b28SMark { 59762c23b28SMark PetscFunctionBegin; 59862c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 59962c23b28SMark if (scale) PetscValidRealPointer(scale,2); 60062c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 60162c23b28SMark PetscFunctionReturn(0); 60262c23b28SMark } 60362c23b28SMark 60462c23b28SMark /*@ 6051917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 6061917a363SLisandro Dalcin 6071917a363SLisandro Dalcin Logically collective on TSAdapt 6081c3436cfSJed Brown 6094165533cSJose E. Roman Input Parameters: 610552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6111c3436cfSJed Brown . hmin - minimum time step 6121c3436cfSJed Brown - hmax - maximum time step 6131c3436cfSJed Brown 6141c3436cfSJed Brown Options Database Keys: 615ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 616ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 6171c3436cfSJed Brown 6181c3436cfSJed Brown Level: intermediate 6191c3436cfSJed Brown 620db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 6211c3436cfSJed Brown @*/ 6221c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 6231c3436cfSJed Brown { 6241c3436cfSJed Brown PetscFunctionBegin; 6254782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 626b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 627b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 6283c633725SBarry Smith PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 6293c633725SBarry Smith PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 630b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 631b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 632b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 633b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6343c633725SBarry 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); 6351917a363SLisandro Dalcin PetscFunctionReturn(0); 6361917a363SLisandro Dalcin } 6371917a363SLisandro Dalcin 6381917a363SLisandro Dalcin /*@ 6391917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 6401917a363SLisandro Dalcin 6411917a363SLisandro Dalcin Not Collective 6421917a363SLisandro Dalcin 6434165533cSJose E. Roman Input Parameter: 6441917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6451917a363SLisandro Dalcin 6464165533cSJose E. Roman Output Parameters: 6471917a363SLisandro Dalcin + hmin - minimum time step 6481917a363SLisandro Dalcin - hmax - maximum time step 6491917a363SLisandro Dalcin 6501917a363SLisandro Dalcin Level: intermediate 6511917a363SLisandro Dalcin 652db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 6531917a363SLisandro Dalcin @*/ 6541917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 6551917a363SLisandro Dalcin { 6561917a363SLisandro Dalcin PetscFunctionBegin; 6571917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6581917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 6591917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 6601917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6611917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6621c3436cfSJed Brown PetscFunctionReturn(0); 6631c3436cfSJed Brown } 6641c3436cfSJed Brown 665e55864a3SBarry Smith /* 66684df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 66784df9cb4SJed Brown 66884df9cb4SJed Brown Collective on TSAdapt 66984df9cb4SJed Brown 67084df9cb4SJed Brown Input Parameter: 67184df9cb4SJed Brown . adapt - the TSAdapt context 67284df9cb4SJed Brown 67384df9cb4SJed Brown Options Database Keys: 6741917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 675bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6761917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6771917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6781917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 67923de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 68023de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 68123de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 682de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 683de50f1caSBarry 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 68484df9cb4SJed Brown 68584df9cb4SJed Brown Level: advanced 68684df9cb4SJed Brown 68784df9cb4SJed Brown Notes: 68884df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 68984df9cb4SJed Brown 690db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 691db781477SPatrick Sanan `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 692e55864a3SBarry Smith */ 6934416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 69484df9cb4SJed Brown { 69584df9cb4SJed Brown char type[256] = TSADAPTBASIC; 69662c23b28SMark PetscReal safety,reject_safety,clip[2],scale,hmin,hmax; 6971c3436cfSJed Brown PetscBool set,flg; 6981917a363SLisandro Dalcin PetscInt two; 69984df9cb4SJed Brown 70084df9cb4SJed Brown PetscFunctionBegin; 701064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,2); 70284df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 7031566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 704d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"TS Adaptivity options"); 7059566063dSJacob 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)); 70684df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 7079566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt,type)); 70884df9cb4SJed Brown } 7091917a363SLisandro Dalcin 7109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set)); 7119566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt,flg)); 7121917a363SLisandro Dalcin 7131917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 7149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set)); 7159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg)); 7169566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt,safety,reject_safety)); 7171917a363SLisandro Dalcin 7181917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 7199566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set)); 7203c633725SBarry Smith PetscCheck(!set || (two == 2),PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 7219566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt,clip[0],clip[1])); 7221917a363SLisandro Dalcin 7231917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 7249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set)); 7259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg)); 7269566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt,hmin,hmax)); 7271917a363SLisandro Dalcin 7289566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set)); 7299566063dSJacob 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)); 730d580f011SEmil Constantinescu 7319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","TSAdaptSetScaleSolveFailed",adapt->scale_solve_failed,&scale,&set)); 7329566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt,scale)); 7331917a363SLisandro Dalcin 7349566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL)); 7353c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY,PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 7361917a363SLisandro Dalcin 7379566063dSJacob 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)); 738de50f1caSBarry Smith 7399566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set)); 7409566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt,flg)); 7411917a363SLisandro Dalcin 7429566063dSJacob Faibussowitsch if (adapt->ops->setfromoptions) PetscCall((*adapt->ops->setfromoptions)(PetscOptionsObject,adapt)); 743d0609cedSBarry Smith PetscOptionsHeadEnd(); 74484df9cb4SJed Brown PetscFunctionReturn(0); 74584df9cb4SJed Brown } 74684df9cb4SJed Brown 74784df9cb4SJed Brown /*@ 74884df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 74984df9cb4SJed Brown 7501917a363SLisandro Dalcin Logically collective on TSAdapt 75184df9cb4SJed Brown 7524165533cSJose E. Roman Input Parameter: 75384df9cb4SJed Brown . adapt - adaptive controller 75484df9cb4SJed Brown 75584df9cb4SJed Brown Level: developer 75684df9cb4SJed Brown 757db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 75884df9cb4SJed Brown @*/ 75984df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 76084df9cb4SJed Brown { 76184df9cb4SJed Brown PetscFunctionBegin; 7624782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7639566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates,sizeof(adapt->candidates))); 76484df9cb4SJed Brown PetscFunctionReturn(0); 76584df9cb4SJed Brown } 76684df9cb4SJed Brown 76784df9cb4SJed Brown /*@C 76884df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 76984df9cb4SJed Brown 7701917a363SLisandro Dalcin Logically collective on TSAdapt 77184df9cb4SJed Brown 7724165533cSJose E. Roman Input Parameters: 773552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 77484df9cb4SJed Brown . name - name of the candidate scheme to add 77584df9cb4SJed Brown . order - order of the candidate scheme 77684df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7778d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 77884df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 77984df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 78084df9cb4SJed Brown 78184df9cb4SJed Brown Note: 78284df9cb4SJed Brown This routine is not available in Fortran. 78384df9cb4SJed Brown 78484df9cb4SJed Brown Level: developer 78584df9cb4SJed Brown 786db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 78784df9cb4SJed Brown @*/ 7888d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 78984df9cb4SJed Brown { 79084df9cb4SJed Brown PetscInt c; 79184df9cb4SJed Brown 79284df9cb4SJed Brown PetscFunctionBegin; 79384df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 79463a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %" PetscInt_FMT " must be a positive integer",order); 79584df9cb4SJed Brown if (inuse) { 7963c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 79784df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 79884df9cb4SJed Brown } 7991c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 8001c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 801bbd56ea5SKarl Rupp 80284df9cb4SJed Brown adapt->candidates.name[c] = name; 80384df9cb4SJed Brown adapt->candidates.order[c] = order; 80484df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 8058d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 80684df9cb4SJed Brown adapt->candidates.cost[c] = cost; 80784df9cb4SJed Brown adapt->candidates.n++; 80884df9cb4SJed Brown PetscFunctionReturn(0); 80984df9cb4SJed Brown } 81084df9cb4SJed Brown 8118d59e960SJed Brown /*@C 8128d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8138d59e960SJed Brown 8148d59e960SJed Brown Not Collective 8158d59e960SJed Brown 8164165533cSJose E. Roman Input Parameter: 8178d59e960SJed Brown . adapt - time step adaptivity context 8188d59e960SJed Brown 8194165533cSJose E. Roman Output Parameters: 8208d59e960SJed Brown + n - number of candidate schemes, always at least 1 8218d59e960SJed Brown . order - the order of each candidate scheme 8228d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8238d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8248d59e960SJed Brown - cost - the relative cost of each scheme 8258d59e960SJed Brown 8268d59e960SJed Brown Level: developer 8278d59e960SJed Brown 8288d59e960SJed Brown Note: 8298d59e960SJed Brown The current scheme is always returned in the first slot 8308d59e960SJed Brown 831db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 8328d59e960SJed Brown @*/ 8338d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 8348d59e960SJed Brown { 8358d59e960SJed Brown PetscFunctionBegin; 8368d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8378d59e960SJed Brown if (n) *n = adapt->candidates.n; 8388d59e960SJed Brown if (order) *order = adapt->candidates.order; 8398d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8408d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8418d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8428d59e960SJed Brown PetscFunctionReturn(0); 8438d59e960SJed Brown } 8448d59e960SJed Brown 84584df9cb4SJed Brown /*@C 84684df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 84784df9cb4SJed Brown 8481917a363SLisandro Dalcin Collective on TSAdapt 84984df9cb4SJed Brown 8504165533cSJose E. Roman Input Parameters: 85184df9cb4SJed Brown + adapt - adaptive contoller 85297bb3fdcSJose E. Roman . ts - time stepper 85384df9cb4SJed Brown - h - current step size 85484df9cb4SJed Brown 8554165533cSJose E. Roman Output Parameters: 8561566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 85784df9cb4SJed Brown . next_h - step size to use for the next step 85884df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 85984df9cb4SJed Brown 8601c3436cfSJed Brown Note: 8611c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 8621c3436cfSJed Brown being retried after an initial rejection. 8631c3436cfSJed Brown 86484df9cb4SJed Brown Level: developer 86584df9cb4SJed Brown 866db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 86784df9cb4SJed Brown @*/ 86884df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 86984df9cb4SJed Brown { 8701566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8711566a47fSLisandro Dalcin PetscInt scheme = 0; 8720b99f514SJed Brown PetscReal wlte = -1.0; 8735e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8745e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 87584df9cb4SJed Brown 87684df9cb4SJed Brown PetscFunctionBegin; 87784df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 87884df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8791566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 880dadcf809SJacob Faibussowitsch PetscValidRealPointer(next_h,5); 881064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,6); 8821566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8831566a47fSLisandro Dalcin 8841566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8851566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8861566a47fSLisandro Dalcin *next_h = h; 8871566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8881566a47fSLisandro Dalcin PetscFunctionReturn(0); 8891566a47fSLisandro Dalcin } 8901566a47fSLisandro Dalcin 8919566063dSJacob Faibussowitsch PetscCall((*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter)); 89263a3b9bcSJacob 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); 8933c633725SBarry Smith PetscCheck(*next_h >= 0,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8941566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8951566a47fSLisandro Dalcin 89649354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 89736b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 89836b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 8994a658b32SHong Zhang PetscReal tend = t + h, tmax, hmax; 900fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 901fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 9024a658b32SHong Zhang 9034a658b32SHong Zhang if (ts->tspan) { 904*e1db57b0SHong Zhang if (PetscIsCloseAtTol(t,ts->tspan->span_times[ts->tspan->spanctr],ts->tspan->reltol*h+ts->tspan->abstol,0)) /* hit a span time point */ 9054a658b32SHong Zhang if (ts->tspan->spanctr+1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr+1]; 9064a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */ 9074a658b32SHong Zhang else tmax = ts->tspan->span_times[ts->tspan->spanctr]; 9084a658b32SHong Zhang } else tmax = ts->max_time; 9094a658b32SHong Zhang hmax = tmax - t; 91036b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 911fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 91236b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 913*e1db57b0SHong Zhang /* if step size is changed to match a span time point */ 914*e1db57b0SHong Zhang if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h; 915*e1db57b0SHong Zhang /* reset time step after a span time point */ 916*e1db57b0SHong 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)) { 917*e1db57b0SHong Zhang *next_h = adapt->dt_span_cached; 918*e1db57b0SHong Zhang adapt->dt_span_cached = 0; 91949354f04SShri Abhyankar } 920*e1db57b0SHong Zhang } 9211c3436cfSJed Brown if (adapt->monitor) { 9221566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 9239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 9240b99f514SJed Brown if (wlte < 0) { 92563a3b9bcSJacob Faibussowitsch 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",(double)ts->ptime,(double)h,(double)*next_h)); 9260b99f514SJed Brown } else { 92763a3b9bcSJacob Faibussowitsch 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",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter)); 9280b99f514SJed Brown } 9299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 9301c3436cfSJed Brown } 93184df9cb4SJed Brown PetscFunctionReturn(0); 93284df9cb4SJed Brown } 93384df9cb4SJed Brown 93497335746SJed Brown /*@ 935de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 936de50f1caSBarry Smith before increasing the time step. 937de50f1caSBarry Smith 938de50f1caSBarry Smith Logicially Collective on TSAdapt 939de50f1caSBarry Smith 9404165533cSJose E. Roman Input Parameters: 941de50f1caSBarry Smith + adapt - adaptive controller context 942de50f1caSBarry Smith - cnt - the number of timesteps 943de50f1caSBarry Smith 944de50f1caSBarry Smith Options Database Key: 945de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 946de50f1caSBarry Smith 947de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 948de50f1caSBarry Smith The successful use of this option is problem dependent 949de50f1caSBarry Smith 950de50f1caSBarry Smith Developer Note: there is no theory to support this option 951de50f1caSBarry Smith 952de50f1caSBarry Smith Level: advanced 953de50f1caSBarry Smith 954de50f1caSBarry Smith .seealso: 955de50f1caSBarry Smith @*/ 956de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 957de50f1caSBarry Smith { 958de50f1caSBarry Smith PetscFunctionBegin; 959de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 960de50f1caSBarry Smith PetscFunctionReturn(0); 961de50f1caSBarry Smith } 962de50f1caSBarry Smith 963de50f1caSBarry Smith /*@ 9646bc98fa9SBarry 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) 96597335746SJed Brown 9661917a363SLisandro Dalcin Collective on TSAdapt 96797335746SJed Brown 9684165533cSJose E. Roman Input Parameters: 96997335746SJed Brown + adapt - adaptive controller context 970b295832fSPierre Barbier de Reuille . ts - time stepper 971b295832fSPierre Barbier de Reuille . t - Current simulation time 972b295832fSPierre Barbier de Reuille - Y - Current solution vector 97397335746SJed Brown 9744165533cSJose E. Roman Output Parameter: 97597335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 97697335746SJed Brown 97797335746SJed Brown Level: developer 97897335746SJed Brown 97997335746SJed Brown .seealso: 98097335746SJed Brown @*/ 981b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 98297335746SJed Brown { 9831566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 98497335746SJed Brown 98597335746SJed Brown PetscFunctionBegin; 9864782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 9874782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 988064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,5); 9891566a47fSLisandro Dalcin 9909566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes,&snesreason)); 99197335746SJed Brown if (snesreason < 0) { 99297335746SJed Brown *accept = PETSC_FALSE; 9936de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 99497335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 99563a3b9bcSJacob 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)); 99697335746SJed Brown if (adapt->monitor) { 9979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 99863a3b9bcSJacob Faibussowitsch 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,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures)); 9999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 100097335746SJed Brown } 1001cb9d8021SPierre Barbier de Reuille } 1002cb9d8021SPierre Barbier de Reuille } else { 10031566a47fSLisandro Dalcin *accept = PETSC_TRUE; 10049566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts,t,Y,accept)); 1005cb9d8021SPierre Barbier de Reuille if (*accept && adapt->checkstage) { 10069566063dSJacob Faibussowitsch PetscCall((*adapt->checkstage)(adapt,ts,t,Y,accept)); 10076bc98fa9SBarry Smith if (!*accept) { 100863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n",ts->steps)); 10096bc98fa9SBarry Smith if (adapt->monitor) { 10109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 101163a3b9bcSJacob 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)); 10129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 10136bc98fa9SBarry Smith } 10146bc98fa9SBarry Smith } 1015cb9d8021SPierre Barbier de Reuille } 1016cb9d8021SPierre Barbier de Reuille } 1017cb9d8021SPierre Barbier de Reuille 10181566a47fSLisandro Dalcin if (!(*accept) && !ts->reason) { 10191566a47fSLisandro Dalcin PetscReal dt,new_dt; 10209566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 1021cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 10229566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,new_dt)); 1023de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 102497335746SJed Brown if (adapt->monitor) { 10259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 102663a3b9bcSJacob 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)); 10279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 102897335746SJed Brown } 102997335746SJed Brown } 103097335746SJed Brown PetscFunctionReturn(0); 103197335746SJed Brown } 103297335746SJed Brown 103384df9cb4SJed Brown /*@ 103484df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 103584df9cb4SJed Brown 1036d083f849SBarry Smith Collective 103784df9cb4SJed Brown 103884df9cb4SJed Brown Input Parameter: 103984df9cb4SJed Brown . comm - The communicator 104084df9cb4SJed Brown 104184df9cb4SJed Brown Output Parameter: 104284df9cb4SJed Brown . adapt - new TSAdapt object 104384df9cb4SJed Brown 104484df9cb4SJed Brown Level: developer 104584df9cb4SJed Brown 104684df9cb4SJed Brown Notes: 104784df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 104884df9cb4SJed Brown 1049db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 105084df9cb4SJed Brown @*/ 105184df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 105284df9cb4SJed Brown { 105384df9cb4SJed Brown TSAdapt adapt; 105484df9cb4SJed Brown 105584df9cb4SJed Brown PetscFunctionBegin; 1056064a246eSJacob Faibussowitsch PetscValidPointer(inadapt,2); 10573b3bcf4cSLisandro Dalcin *inadapt = NULL; 10589566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 10593b3bcf4cSLisandro Dalcin 10609566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView)); 10611c3436cfSJed Brown 1062bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 10631917a363SLisandro Dalcin adapt->safety = 0.9; 10641917a363SLisandro Dalcin adapt->reject_safety = 0.5; 10651917a363SLisandro Dalcin adapt->clip[0] = 0.1; 10661917a363SLisandro Dalcin adapt->clip[1] = 10.; 10671c3436cfSJed Brown adapt->dt_min = 1e-20; 10681917a363SLisandro Dalcin adapt->dt_max = 1e+20; 10691c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1070d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 107197335746SJed Brown adapt->scale_solve_failed = 0.25; 1072fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1073fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1074fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1075fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 10767619abb3SShri adapt->wnormtype = NORM_2; 1077de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 10781917a363SLisandro Dalcin 107984df9cb4SJed Brown *inadapt = adapt; 108084df9cb4SJed Brown PetscFunctionReturn(0); 108184df9cb4SJed Brown } 1082