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 411c84c290SBarry Smith .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 5884df9cb4SJed Brown .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 8084df9cb4SJed Brown .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 9784df9cb4SJed Brown .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 1247eef6055SBarry Smith .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 158d0288e4fSLisandro Dalcin .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 202ad6bc421SBarry Smith .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)); 264bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 2659566063dSJacob Faibussowitsch PetscCall((*adapt->ops->view)(adapt,viewer)); 266f2c2a1b9SBarry Smith } 26784df9cb4SJed Brown PetscFunctionReturn(0); 26884df9cb4SJed Brown } 26984df9cb4SJed Brown 270099af0a3SLisandro Dalcin /*@ 271099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 272099af0a3SLisandro Dalcin 273099af0a3SLisandro Dalcin Collective on TS 274099af0a3SLisandro Dalcin 275099af0a3SLisandro Dalcin Input Parameter: 276099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 277099af0a3SLisandro Dalcin 278099af0a3SLisandro Dalcin Level: developer 279099af0a3SLisandro Dalcin 280099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 281099af0a3SLisandro Dalcin @*/ 282099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 283099af0a3SLisandro Dalcin { 284099af0a3SLisandro Dalcin PetscFunctionBegin; 285099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2869566063dSJacob Faibussowitsch if (adapt->ops->reset) PetscCall((*adapt->ops->reset)(adapt)); 287099af0a3SLisandro Dalcin PetscFunctionReturn(0); 288099af0a3SLisandro Dalcin } 289099af0a3SLisandro Dalcin 29084df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 29184df9cb4SJed Brown { 29284df9cb4SJed Brown PetscFunctionBegin; 29384df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 29484df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 2954782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 296099af0a3SLisandro Dalcin 2979566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt)); 298099af0a3SLisandro Dalcin 2999566063dSJacob Faibussowitsch if ((*adapt)->ops->destroy) PetscCall((*(*adapt)->ops->destroy)(*adapt)); 3009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 3019566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt)); 30284df9cb4SJed Brown PetscFunctionReturn(0); 30384df9cb4SJed Brown } 30484df9cb4SJed Brown 3051c3436cfSJed Brown /*@ 3061c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3071c3436cfSJed Brown 3081c3436cfSJed Brown Collective on TSAdapt 3091c3436cfSJed Brown 3104165533cSJose E. Roman Input Parameters: 3111c3436cfSJed Brown + adapt - adaptive controller context 3121c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3131c3436cfSJed Brown 314bf997491SLisandro Dalcin Options Database Keys: 315ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 316bf997491SLisandro Dalcin 3171c3436cfSJed Brown Level: intermediate 3181c3436cfSJed Brown 3191c3436cfSJed Brown .seealso: TSAdaptChoose() 3201c3436cfSJed Brown @*/ 3211c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3221c3436cfSJed Brown { 3231c3436cfSJed Brown PetscFunctionBegin; 3244782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3254782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3261c3436cfSJed Brown if (flg) { 3279566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor)); 3281c3436cfSJed Brown } else { 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor)); 3301c3436cfSJed Brown } 3311c3436cfSJed Brown PetscFunctionReturn(0); 3321c3436cfSJed Brown } 3331c3436cfSJed Brown 3340873808bSJed Brown /*@C 335bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3360873808bSJed Brown 3370873808bSJed Brown Logically collective on TSAdapt 3380873808bSJed Brown 3394165533cSJose E. Roman Input Parameters: 3400873808bSJed Brown + adapt - adaptive controller context 3410873808bSJed Brown - func - stage check function 3420873808bSJed Brown 3430873808bSJed Brown Arguments of func: 3440873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3450873808bSJed Brown 3460873808bSJed Brown + adapt - adaptive controller context 3470873808bSJed Brown . ts - time stepping context 3480873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3490873808bSJed Brown 3500873808bSJed Brown Level: advanced 3510873808bSJed Brown 3520873808bSJed Brown .seealso: TSAdaptChoose() 3530873808bSJed Brown @*/ 354b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3550873808bSJed Brown { 3560873808bSJed Brown PetscFunctionBegin; 3570873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 35868ae941aSLisandro Dalcin adapt->checkstage = func; 3590873808bSJed Brown PetscFunctionReturn(0); 3600873808bSJed Brown } 3610873808bSJed Brown 3621c3436cfSJed Brown /*@ 363bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 364bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 365bf997491SLisandro Dalcin 3661917a363SLisandro Dalcin Logically collective on TSAdapt 367bf997491SLisandro Dalcin 3684165533cSJose E. Roman Input Parameters: 369bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 370bf997491SLisandro Dalcin - flag - whether to always accept steps 371bf997491SLisandro Dalcin 372bf997491SLisandro Dalcin Options Database Keys: 373ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 374bf997491SLisandro Dalcin 375bf997491SLisandro Dalcin Level: intermediate 376bf997491SLisandro Dalcin 377bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 378bf997491SLisandro Dalcin @*/ 379bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 380bf997491SLisandro Dalcin { 381bf997491SLisandro Dalcin PetscFunctionBegin; 382bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 383bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 384bf997491SLisandro Dalcin adapt->always_accept = flag; 385bf997491SLisandro Dalcin PetscFunctionReturn(0); 386bf997491SLisandro Dalcin } 387bf997491SLisandro Dalcin 388bf997491SLisandro Dalcin /*@ 3891917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 3901c3436cfSJed Brown 3911917a363SLisandro Dalcin Logically collective on TSAdapt 3921917a363SLisandro Dalcin 3934165533cSJose E. Roman Input Parameters: 3941917a363SLisandro Dalcin + adapt - adaptive controller context 3951917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 396ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 3971917a363SLisandro Dalcin 3981917a363SLisandro Dalcin Options Database Keys: 399ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 400ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4011917a363SLisandro Dalcin 4021917a363SLisandro Dalcin Level: intermediate 4031917a363SLisandro Dalcin 4041917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4051917a363SLisandro Dalcin @*/ 4061917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4071917a363SLisandro Dalcin { 4081917a363SLisandro Dalcin PetscFunctionBegin; 4091917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4101917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4111917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4123c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4133c633725SBarry Smith PetscCheck(safety == PETSC_DEFAULT || safety <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety); 4143c633725SBarry 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); 4153c633725SBarry 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); 4161917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4171917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4181917a363SLisandro Dalcin PetscFunctionReturn(0); 4191917a363SLisandro Dalcin } 4201917a363SLisandro Dalcin 4211917a363SLisandro Dalcin /*@ 4221917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4231917a363SLisandro Dalcin 4241917a363SLisandro Dalcin Not Collective 4251917a363SLisandro Dalcin 4264165533cSJose E. Roman Input Parameter: 4271917a363SLisandro Dalcin . adapt - adaptive controller context 4281917a363SLisandro Dalcin 4294165533cSJose E. Roman Output Parameters: 4301917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4311917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4321917a363SLisandro Dalcin 4331917a363SLisandro Dalcin Level: intermediate 4341917a363SLisandro Dalcin 4351917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4361917a363SLisandro Dalcin @*/ 4371917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4381917a363SLisandro Dalcin { 4391917a363SLisandro Dalcin PetscFunctionBegin; 4401917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4411917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4421917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4431917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4441917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4451917a363SLisandro Dalcin PetscFunctionReturn(0); 4461917a363SLisandro Dalcin } 4471917a363SLisandro Dalcin 4481917a363SLisandro Dalcin /*@ 44976cddca1SEmil 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. 45076cddca1SEmil Constantinescu 45176cddca1SEmil Constantinescu Logically collective on TSAdapt 45276cddca1SEmil Constantinescu 4534165533cSJose E. Roman Input Parameters: 45476cddca1SEmil Constantinescu + adapt - adaptive controller context 45576cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 45676cddca1SEmil Constantinescu 45776cddca1SEmil Constantinescu Options Database Keys: 45876cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 45976cddca1SEmil Constantinescu 46076cddca1SEmil Constantinescu Level: intermediate 46176cddca1SEmil Constantinescu 46276cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptGetMaxIgnore(), TSAdaptChoose() 46376cddca1SEmil Constantinescu @*/ 46476cddca1SEmil Constantinescu PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore) 46576cddca1SEmil Constantinescu { 46676cddca1SEmil Constantinescu PetscFunctionBegin; 46776cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 46876cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt,max_ignore,2); 46976cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 47076cddca1SEmil Constantinescu PetscFunctionReturn(0); 47176cddca1SEmil Constantinescu } 47276cddca1SEmil Constantinescu 47376cddca1SEmil Constantinescu /*@ 47476cddca1SEmil 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). 47576cddca1SEmil Constantinescu 47676cddca1SEmil Constantinescu Not Collective 47776cddca1SEmil Constantinescu 4784165533cSJose E. Roman Input Parameter: 47976cddca1SEmil Constantinescu . adapt - adaptive controller context 48076cddca1SEmil Constantinescu 4814165533cSJose E. Roman Output Parameter: 48276cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 48376cddca1SEmil Constantinescu 48476cddca1SEmil Constantinescu Level: intermediate 48576cddca1SEmil Constantinescu 48676cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptSetMaxIgnore(), TSAdaptChoose() 48776cddca1SEmil Constantinescu @*/ 48876cddca1SEmil Constantinescu PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore) 48976cddca1SEmil Constantinescu { 49076cddca1SEmil Constantinescu PetscFunctionBegin; 49176cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 49276cddca1SEmil Constantinescu PetscValidRealPointer(max_ignore,2); 49376cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 49476cddca1SEmil Constantinescu PetscFunctionReturn(0); 49576cddca1SEmil Constantinescu } 49676cddca1SEmil Constantinescu 49776cddca1SEmil Constantinescu /*@ 4981917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4991917a363SLisandro Dalcin 5001917a363SLisandro Dalcin Logically collective on TSAdapt 5011917a363SLisandro Dalcin 5024165533cSJose E. Roman Input Parameters: 5031917a363SLisandro Dalcin + adapt - adaptive controller context 5041917a363SLisandro Dalcin . low - admissible decrease factor 505ec18b7bcSLisandro Dalcin - high - admissible increase factor 5061917a363SLisandro Dalcin 5071917a363SLisandro Dalcin Options Database Keys: 508ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5091917a363SLisandro Dalcin 5101917a363SLisandro Dalcin Level: intermediate 5111917a363SLisandro Dalcin 51262c23b28SMark .seealso: TSAdaptChoose(), TSAdaptGetClip(), TSAdaptSetScaleSolveFailed() 5131917a363SLisandro Dalcin @*/ 5141917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 5151917a363SLisandro Dalcin { 5161917a363SLisandro Dalcin PetscFunctionBegin; 5171917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5181917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 5191917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 5203c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 5213c633725SBarry Smith PetscCheck(low == PETSC_DEFAULT || low <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 5223c633725SBarry Smith PetscCheck(high == PETSC_DEFAULT || high >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be greater than one",(double)high); 5231917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 5241917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 5251917a363SLisandro Dalcin PetscFunctionReturn(0); 5261917a363SLisandro Dalcin } 5271917a363SLisandro Dalcin 5281917a363SLisandro Dalcin /*@ 5291917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5301917a363SLisandro Dalcin 5311917a363SLisandro Dalcin Not Collective 5321917a363SLisandro Dalcin 5334165533cSJose E. Roman Input Parameter: 5341917a363SLisandro Dalcin . adapt - adaptive controller context 5351917a363SLisandro Dalcin 5364165533cSJose E. Roman Output Parameters: 5371917a363SLisandro Dalcin + low - optional, admissible decrease factor 5381917a363SLisandro Dalcin - high - optional, admissible increase factor 5391917a363SLisandro Dalcin 5401917a363SLisandro Dalcin Level: intermediate 5411917a363SLisandro Dalcin 54262c23b28SMark .seealso: TSAdaptChoose(), TSAdaptSetClip(), TSAdaptSetScaleSolveFailed() 5431917a363SLisandro Dalcin @*/ 5441917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5451917a363SLisandro Dalcin { 5461917a363SLisandro Dalcin PetscFunctionBegin; 5471917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5481917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5491917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5501917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5511917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5521917a363SLisandro Dalcin PetscFunctionReturn(0); 5531917a363SLisandro Dalcin } 5541917a363SLisandro Dalcin 5551917a363SLisandro Dalcin /*@ 55662c23b28SMark TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 55762c23b28SMark 55862c23b28SMark Logically collective on TSAdapt 55962c23b28SMark 5604165533cSJose E. Roman Input Parameters: 56162c23b28SMark + adapt - adaptive controller context 562ee300463SSatish Balay - scale - scale 56362c23b28SMark 56462c23b28SMark Options Database Keys: 56562c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 56662c23b28SMark 56762c23b28SMark Level: intermediate 56862c23b28SMark 56962c23b28SMark .seealso: TSAdaptChoose(), TSAdaptGetScaleSolveFailed(), TSAdaptGetClip() 57062c23b28SMark @*/ 57162c23b28SMark PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt,PetscReal scale) 57262c23b28SMark { 57362c23b28SMark PetscFunctionBegin; 57462c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 57562c23b28SMark PetscValidLogicalCollectiveReal(adapt,scale,2); 5763c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be positive",(double)scale); 5773c633725SBarry Smith PetscCheck(scale == PETSC_DEFAULT || scale <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be less than one",(double)scale); 57862c23b28SMark if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 57962c23b28SMark PetscFunctionReturn(0); 58062c23b28SMark } 58162c23b28SMark 58262c23b28SMark /*@ 58362c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 58462c23b28SMark 58562c23b28SMark Not Collective 58662c23b28SMark 5874165533cSJose E. Roman Input Parameter: 58862c23b28SMark . adapt - adaptive controller context 58962c23b28SMark 5904165533cSJose E. Roman Output Parameter: 591ee300463SSatish Balay . scale - scale factor 59262c23b28SMark 59362c23b28SMark Level: intermediate 59462c23b28SMark 59562c23b28SMark .seealso: TSAdaptChoose(), TSAdaptSetScaleSolveFailed(), TSAdaptSetClip() 59662c23b28SMark @*/ 59762c23b28SMark PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt,PetscReal *scale) 59862c23b28SMark { 59962c23b28SMark PetscFunctionBegin; 60062c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 60162c23b28SMark if (scale) PetscValidRealPointer(scale,2); 60262c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 60362c23b28SMark PetscFunctionReturn(0); 60462c23b28SMark } 60562c23b28SMark 60662c23b28SMark /*@ 6071917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 6081917a363SLisandro Dalcin 6091917a363SLisandro Dalcin Logically collective on TSAdapt 6101c3436cfSJed Brown 6114165533cSJose E. Roman Input Parameters: 612552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6131c3436cfSJed Brown . hmin - minimum time step 6141c3436cfSJed Brown - hmax - maximum time step 6151c3436cfSJed Brown 6161c3436cfSJed Brown Options Database Keys: 617ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 618ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 6191c3436cfSJed Brown 6201c3436cfSJed Brown Level: intermediate 6211c3436cfSJed Brown 6221917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 6231c3436cfSJed Brown @*/ 6241c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 6251c3436cfSJed Brown { 6261c3436cfSJed Brown PetscFunctionBegin; 6274782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 628b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 629b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 6303c633725SBarry Smith PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 6313c633725SBarry Smith PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 632b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 633b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 634b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 635b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6363c633725SBarry 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); 6371917a363SLisandro Dalcin PetscFunctionReturn(0); 6381917a363SLisandro Dalcin } 6391917a363SLisandro Dalcin 6401917a363SLisandro Dalcin /*@ 6411917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 6421917a363SLisandro Dalcin 6431917a363SLisandro Dalcin Not Collective 6441917a363SLisandro Dalcin 6454165533cSJose E. Roman Input Parameter: 6461917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6471917a363SLisandro Dalcin 6484165533cSJose E. Roman Output Parameters: 6491917a363SLisandro Dalcin + hmin - minimum time step 6501917a363SLisandro Dalcin - hmax - maximum time step 6511917a363SLisandro Dalcin 6521917a363SLisandro Dalcin Level: intermediate 6531917a363SLisandro Dalcin 6541917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 6551917a363SLisandro Dalcin @*/ 6561917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 6571917a363SLisandro Dalcin { 6581917a363SLisandro Dalcin PetscFunctionBegin; 6591917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6601917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 6611917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 6621917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6631917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6641c3436cfSJed Brown PetscFunctionReturn(0); 6651c3436cfSJed Brown } 6661c3436cfSJed Brown 667e55864a3SBarry Smith /* 66884df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 66984df9cb4SJed Brown 67084df9cb4SJed Brown Collective on TSAdapt 67184df9cb4SJed Brown 67284df9cb4SJed Brown Input Parameter: 67384df9cb4SJed Brown . adapt - the TSAdapt context 67484df9cb4SJed Brown 67584df9cb4SJed Brown Options Database Keys: 6761917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 677bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6781917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6791917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6801917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 68123de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 68223de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 68323de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 684de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 685de50f1caSBarry 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 68684df9cb4SJed Brown 68784df9cb4SJed Brown Level: advanced 68884df9cb4SJed Brown 68984df9cb4SJed Brown Notes: 69084df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 69184df9cb4SJed Brown 6921917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 69362c23b28SMark TSAdaptSetClip(), TSAdaptSetScaleSolveFailed(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 694e55864a3SBarry Smith */ 6954416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 69684df9cb4SJed Brown { 69784df9cb4SJed Brown char type[256] = TSADAPTBASIC; 69862c23b28SMark PetscReal safety,reject_safety,clip[2],scale,hmin,hmax; 6991c3436cfSJed Brown PetscBool set,flg; 7001917a363SLisandro Dalcin PetscInt two; 70184df9cb4SJed Brown 70284df9cb4SJed Brown PetscFunctionBegin; 703064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,2); 70484df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 7051566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 7069566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options")); 7079566063dSJacob 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)); 70884df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 7099566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt,type)); 71084df9cb4SJed Brown } 7111917a363SLisandro Dalcin 7129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set)); 7139566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt,flg)); 7141917a363SLisandro Dalcin 7151917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 7169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set)); 7179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg)); 7189566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt,safety,reject_safety)); 7191917a363SLisandro Dalcin 7201917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 7219566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set)); 7223c633725SBarry Smith PetscCheck(!set || (two == 2),PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 7239566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt,clip[0],clip[1])); 7241917a363SLisandro Dalcin 7251917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 7269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set)); 7279566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg)); 7289566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt,hmin,hmax)); 7291917a363SLisandro Dalcin 7309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set)); 7319566063dSJacob 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)); 732d580f011SEmil Constantinescu 7339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","TSAdaptSetScaleSolveFailed",adapt->scale_solve_failed,&scale,&set)); 7349566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt,scale)); 7351917a363SLisandro Dalcin 7369566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL)); 7373c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY,PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 7381917a363SLisandro Dalcin 7399566063dSJacob 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)); 740de50f1caSBarry Smith 7419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set)); 7429566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt,flg)); 7431917a363SLisandro Dalcin 7449566063dSJacob Faibussowitsch if (adapt->ops->setfromoptions) PetscCall((*adapt->ops->setfromoptions)(PetscOptionsObject,adapt)); 7459566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 74684df9cb4SJed Brown PetscFunctionReturn(0); 74784df9cb4SJed Brown } 74884df9cb4SJed Brown 74984df9cb4SJed Brown /*@ 75084df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 75184df9cb4SJed Brown 7521917a363SLisandro Dalcin Logically collective on TSAdapt 75384df9cb4SJed Brown 7544165533cSJose E. Roman Input Parameter: 75584df9cb4SJed Brown . adapt - adaptive controller 75684df9cb4SJed Brown 75784df9cb4SJed Brown Level: developer 75884df9cb4SJed Brown 75984df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 76084df9cb4SJed Brown @*/ 76184df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 76284df9cb4SJed Brown { 76384df9cb4SJed Brown PetscFunctionBegin; 7644782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7659566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates,sizeof(adapt->candidates))); 76684df9cb4SJed Brown PetscFunctionReturn(0); 76784df9cb4SJed Brown } 76884df9cb4SJed Brown 76984df9cb4SJed Brown /*@C 77084df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 77184df9cb4SJed Brown 7721917a363SLisandro Dalcin Logically collective on TSAdapt 77384df9cb4SJed Brown 7744165533cSJose E. Roman Input Parameters: 775552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 77684df9cb4SJed Brown . name - name of the candidate scheme to add 77784df9cb4SJed Brown . order - order of the candidate scheme 77884df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7798d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 78084df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 78184df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 78284df9cb4SJed Brown 78384df9cb4SJed Brown Note: 78484df9cb4SJed Brown This routine is not available in Fortran. 78584df9cb4SJed Brown 78684df9cb4SJed Brown Level: developer 78784df9cb4SJed Brown 78884df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 78984df9cb4SJed Brown @*/ 7908d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 79184df9cb4SJed Brown { 79284df9cb4SJed Brown PetscInt c; 79384df9cb4SJed Brown 79484df9cb4SJed Brown PetscFunctionBegin; 79584df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7963c633725SBarry Smith PetscCheck(order >= 1,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 79784df9cb4SJed Brown if (inuse) { 7983c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 79984df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 80084df9cb4SJed Brown } 8011c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 8021c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 803bbd56ea5SKarl Rupp 80484df9cb4SJed Brown adapt->candidates.name[c] = name; 80584df9cb4SJed Brown adapt->candidates.order[c] = order; 80684df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 8078d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 80884df9cb4SJed Brown adapt->candidates.cost[c] = cost; 80984df9cb4SJed Brown adapt->candidates.n++; 81084df9cb4SJed Brown PetscFunctionReturn(0); 81184df9cb4SJed Brown } 81284df9cb4SJed Brown 8138d59e960SJed Brown /*@C 8148d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8158d59e960SJed Brown 8168d59e960SJed Brown Not Collective 8178d59e960SJed Brown 8184165533cSJose E. Roman Input Parameter: 8198d59e960SJed Brown . adapt - time step adaptivity context 8208d59e960SJed Brown 8214165533cSJose E. Roman Output Parameters: 8228d59e960SJed Brown + n - number of candidate schemes, always at least 1 8238d59e960SJed Brown . order - the order of each candidate scheme 8248d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8258d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8268d59e960SJed Brown - cost - the relative cost of each scheme 8278d59e960SJed Brown 8288d59e960SJed Brown Level: developer 8298d59e960SJed Brown 8308d59e960SJed Brown Note: 8318d59e960SJed Brown The current scheme is always returned in the first slot 8328d59e960SJed Brown 8338d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 8348d59e960SJed Brown @*/ 8358d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 8368d59e960SJed Brown { 8378d59e960SJed Brown PetscFunctionBegin; 8388d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8398d59e960SJed Brown if (n) *n = adapt->candidates.n; 8408d59e960SJed Brown if (order) *order = adapt->candidates.order; 8418d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8428d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8438d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8448d59e960SJed Brown PetscFunctionReturn(0); 8458d59e960SJed Brown } 8468d59e960SJed Brown 84784df9cb4SJed Brown /*@C 84884df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 84984df9cb4SJed Brown 8501917a363SLisandro Dalcin Collective on TSAdapt 85184df9cb4SJed Brown 8524165533cSJose E. Roman Input Parameters: 85384df9cb4SJed Brown + adapt - adaptive contoller 85497bb3fdcSJose E. Roman . ts - time stepper 85584df9cb4SJed Brown - h - current step size 85684df9cb4SJed Brown 8574165533cSJose E. Roman Output Parameters: 8581566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 85984df9cb4SJed Brown . next_h - step size to use for the next step 86084df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 86184df9cb4SJed Brown 8621c3436cfSJed Brown Note: 8631c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 8641c3436cfSJed Brown being retried after an initial rejection. 8651c3436cfSJed Brown 86684df9cb4SJed Brown Level: developer 86784df9cb4SJed Brown 86884df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 86984df9cb4SJed Brown @*/ 87084df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 87184df9cb4SJed Brown { 8721566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8731566a47fSLisandro Dalcin PetscInt scheme = 0; 8740b99f514SJed Brown PetscReal wlte = -1.0; 8755e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8765e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 87784df9cb4SJed Brown 87884df9cb4SJed Brown PetscFunctionBegin; 87984df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 88084df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8811566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 882dadcf809SJacob Faibussowitsch PetscValidRealPointer(next_h,5); 883064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,6); 8841566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8851566a47fSLisandro Dalcin 8861566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8871566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8881566a47fSLisandro Dalcin *next_h = h; 8891566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8901566a47fSLisandro Dalcin PetscFunctionReturn(0); 8911566a47fSLisandro Dalcin } 8921566a47fSLisandro Dalcin 8939566063dSJacob Faibussowitsch PetscCall((*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter)); 8943c633725SBarry Smith PetscCheck(scheme >= 0 && (ncandidates <= 0 || scheme < ncandidates),PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",scheme,ncandidates-1); 8953c633725SBarry Smith PetscCheck(*next_h >= 0,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8961566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8971566a47fSLisandro Dalcin 89849354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 89936b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 90036b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 901*4a658b32SHong Zhang PetscReal tend = t + h, tmax, hmax; 902fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 903fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 904*4a658b32SHong Zhang 905*4a658b32SHong Zhang if (ts->tspan) { 906*4a658b32SHong Zhang if (PetscIsCloseAtTol(t,ts->tspan->span_times[ts->tspan->spanctr],10*PETSC_MACHINE_EPSILON,0)) /* hit a span time point */ 907*4a658b32SHong Zhang if (ts->tspan->spanctr+1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr+1]; 908*4a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */ 909*4a658b32SHong Zhang else tmax = ts->tspan->span_times[ts->tspan->spanctr]; 910*4a658b32SHong Zhang } else tmax = ts->max_time; 911*4a658b32SHong Zhang hmax = tmax - t; 91236b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 913fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 91436b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 91549354f04SShri Abhyankar } 9161c3436cfSJed Brown 9171c3436cfSJed Brown if (adapt->monitor) { 9181566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 9199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 9200b99f514SJed Brown if (wlte < 0) { 9219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %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)); 9220b99f514SJed Brown } else { 9239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %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)); 9240b99f514SJed Brown } 9259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 9261c3436cfSJed Brown } 92784df9cb4SJed Brown PetscFunctionReturn(0); 92884df9cb4SJed Brown } 92984df9cb4SJed Brown 93097335746SJed Brown /*@ 931de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 932de50f1caSBarry Smith before increasing the time step. 933de50f1caSBarry Smith 934de50f1caSBarry Smith Logicially Collective on TSAdapt 935de50f1caSBarry Smith 9364165533cSJose E. Roman Input Parameters: 937de50f1caSBarry Smith + adapt - adaptive controller context 938de50f1caSBarry Smith - cnt - the number of timesteps 939de50f1caSBarry Smith 940de50f1caSBarry Smith Options Database Key: 941de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 942de50f1caSBarry Smith 943de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 944de50f1caSBarry Smith The successful use of this option is problem dependent 945de50f1caSBarry Smith 946de50f1caSBarry Smith Developer Note: there is no theory to support this option 947de50f1caSBarry Smith 948de50f1caSBarry Smith Level: advanced 949de50f1caSBarry Smith 950de50f1caSBarry Smith .seealso: 951de50f1caSBarry Smith @*/ 952de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 953de50f1caSBarry Smith { 954de50f1caSBarry Smith PetscFunctionBegin; 955de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 956de50f1caSBarry Smith PetscFunctionReturn(0); 957de50f1caSBarry Smith } 958de50f1caSBarry Smith 959de50f1caSBarry Smith /*@ 9606bc98fa9SBarry 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) 96197335746SJed Brown 9621917a363SLisandro Dalcin Collective on TSAdapt 96397335746SJed Brown 9644165533cSJose E. Roman Input Parameters: 96597335746SJed Brown + adapt - adaptive controller context 966b295832fSPierre Barbier de Reuille . ts - time stepper 967b295832fSPierre Barbier de Reuille . t - Current simulation time 968b295832fSPierre Barbier de Reuille - Y - Current solution vector 96997335746SJed Brown 9704165533cSJose E. Roman Output Parameter: 97197335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 97297335746SJed Brown 97397335746SJed Brown Level: developer 97497335746SJed Brown 97597335746SJed Brown .seealso: 97697335746SJed Brown @*/ 977b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 97897335746SJed Brown { 9791566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 98097335746SJed Brown 98197335746SJed Brown PetscFunctionBegin; 9824782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 9834782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 984064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,5); 9851566a47fSLisandro Dalcin 9869566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes,&snesreason)); 98797335746SJed Brown if (snesreason < 0) { 98897335746SJed Brown *accept = PETSC_FALSE; 9896de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 99097335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 9919566063dSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures)); 99297335746SJed Brown if (adapt->monitor) { 9939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 9949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected t=%-11g+%10.3e, nonlinear solve failures %D greater than current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures)); 9959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 99697335746SJed Brown } 997cb9d8021SPierre Barbier de Reuille } 998cb9d8021SPierre Barbier de Reuille } else { 9991566a47fSLisandro Dalcin *accept = PETSC_TRUE; 10009566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts,t,Y,accept)); 1001cb9d8021SPierre Barbier de Reuille if (*accept && adapt->checkstage) { 10029566063dSJacob Faibussowitsch PetscCall((*adapt->checkstage)(adapt,ts,t,Y,accept)); 10036bc98fa9SBarry Smith if (!*accept) { 10049566063dSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%D, solution rejected by user function provided by TSSetFunctionDomainError()\n",ts->steps)); 10056bc98fa9SBarry Smith if (adapt->monitor) { 10069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 10079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected by user function provided by TSSetFunctionDomainError()\n",((PetscObject)adapt)->type_name,ts->steps)); 10089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 10096bc98fa9SBarry Smith } 10106bc98fa9SBarry Smith } 1011cb9d8021SPierre Barbier de Reuille } 1012cb9d8021SPierre Barbier de Reuille } 1013cb9d8021SPierre Barbier de Reuille 10141566a47fSLisandro Dalcin if (!(*accept) && !ts->reason) { 10151566a47fSLisandro Dalcin PetscReal dt,new_dt; 10169566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 1017cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 10189566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,new_dt)); 1019de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 102097335746SJed Brown if (adapt->monitor) { 10219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 10229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D 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)); 10239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 102497335746SJed Brown } 102597335746SJed Brown } 102697335746SJed Brown PetscFunctionReturn(0); 102797335746SJed Brown } 102897335746SJed Brown 102984df9cb4SJed Brown /*@ 103084df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 103184df9cb4SJed Brown 1032d083f849SBarry Smith Collective 103384df9cb4SJed Brown 103484df9cb4SJed Brown Input Parameter: 103584df9cb4SJed Brown . comm - The communicator 103684df9cb4SJed Brown 103784df9cb4SJed Brown Output Parameter: 103884df9cb4SJed Brown . adapt - new TSAdapt object 103984df9cb4SJed Brown 104084df9cb4SJed Brown Level: developer 104184df9cb4SJed Brown 104284df9cb4SJed Brown Notes: 104384df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 104484df9cb4SJed Brown 1045552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 104684df9cb4SJed Brown @*/ 104784df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 104884df9cb4SJed Brown { 104984df9cb4SJed Brown TSAdapt adapt; 105084df9cb4SJed Brown 105184df9cb4SJed Brown PetscFunctionBegin; 1052064a246eSJacob Faibussowitsch PetscValidPointer(inadapt,2); 10533b3bcf4cSLisandro Dalcin *inadapt = NULL; 10549566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage()); 10553b3bcf4cSLisandro Dalcin 10569566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView)); 10571c3436cfSJed Brown 1058bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 10591917a363SLisandro Dalcin adapt->safety = 0.9; 10601917a363SLisandro Dalcin adapt->reject_safety = 0.5; 10611917a363SLisandro Dalcin adapt->clip[0] = 0.1; 10621917a363SLisandro Dalcin adapt->clip[1] = 10.; 10631c3436cfSJed Brown adapt->dt_min = 1e-20; 10641917a363SLisandro Dalcin adapt->dt_max = 1e+20; 10651c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1066d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 106797335746SJed Brown adapt->scale_solve_failed = 0.25; 1068fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1069fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1070fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1071fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 10727619abb3SShri adapt->wnormtype = NORM_2; 1073de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 10741917a363SLisandro Dalcin 107584df9cb4SJed Brown *inadapt = adapt; 107684df9cb4SJed Brown PetscFunctionReturn(0); 107784df9cb4SJed Brown } 1078