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 PetscErrorCode ierr; 4684df9cb4SJed Brown 4784df9cb4SJed Brown PetscFunctionBegin; 481d36bdfdSBarry Smith ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 49a240a19fSJed Brown ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 5084df9cb4SJed Brown PetscFunctionReturn(0); 5184df9cb4SJed Brown } 5284df9cb4SJed Brown 5384df9cb4SJed Brown /*@C 5484df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5584df9cb4SJed Brown 5684df9cb4SJed Brown Not Collective 5784df9cb4SJed Brown 5884df9cb4SJed Brown Level: advanced 5984df9cb4SJed Brown 6084df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6184df9cb4SJed Brown @*/ 62607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6384df9cb4SJed Brown { 6484df9cb4SJed Brown PetscErrorCode ierr; 6584df9cb4SJed Brown 6684df9cb4SJed Brown PetscFunctionBegin; 670f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 680f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 69bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 701566a47fSLisandro Dalcin ierr = TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic);CHKERRQ(ierr); 71cb7ab186SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP);CHKERRQ(ierr); 72bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 731917a363SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE);CHKERRQ(ierr); 74d949e4a4SStefano Zampini ierr = TSAdaptRegister(TSADAPTHISTORY,TSAdaptCreate_History);CHKERRQ(ierr); 7584df9cb4SJed Brown PetscFunctionReturn(0); 7684df9cb4SJed Brown } 7784df9cb4SJed Brown 7884df9cb4SJed Brown /*@C 79560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 8084df9cb4SJed Brown called from PetscFinalize(). 8184df9cb4SJed Brown 8284df9cb4SJed Brown Level: developer 8384df9cb4SJed Brown 8484df9cb4SJed Brown .seealso: PetscFinalize() 8584df9cb4SJed Brown @*/ 8684df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8784df9cb4SJed Brown { 8837e93019SBarry Smith PetscErrorCode ierr; 8937e93019SBarry Smith 9084df9cb4SJed Brown PetscFunctionBegin; 9137e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9284df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9384df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9484df9cb4SJed Brown PetscFunctionReturn(0); 9584df9cb4SJed Brown } 9684df9cb4SJed Brown 9784df9cb4SJed Brown /*@C 9884df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 998a690491SBarry Smith called from TSInitializePackage(). 10084df9cb4SJed Brown 10184df9cb4SJed Brown Level: developer 10284df9cb4SJed Brown 10384df9cb4SJed Brown .seealso: PetscInitialize() 10484df9cb4SJed Brown @*/ 105607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 10684df9cb4SJed Brown { 10784df9cb4SJed Brown PetscErrorCode ierr; 10884df9cb4SJed Brown 10984df9cb4SJed Brown PetscFunctionBegin; 11084df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11184df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11284df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 113607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 11484df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 11584df9cb4SJed Brown PetscFunctionReturn(0); 11684df9cb4SJed Brown } 11784df9cb4SJed Brown 1187eef6055SBarry Smith /*@C 1197eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1207eef6055SBarry Smith 1217eef6055SBarry Smith Logicially Collective on TSAdapt 1227eef6055SBarry Smith 123*d8d19677SJose E. Roman Input Parameters: 124d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1257eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1267eef6055SBarry Smith 1277eef6055SBarry Smith Options Database: 128ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1297eef6055SBarry Smith 1307eef6055SBarry Smith Level: intermediate 1317eef6055SBarry Smith 1327eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1337eef6055SBarry Smith @*/ 13419fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 13584df9cb4SJed Brown { 136ef749922SLisandro Dalcin PetscBool match; 13784df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 13884df9cb4SJed Brown 13984df9cb4SJed Brown PetscFunctionBegin; 1404782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 141b92453a8SLisandro Dalcin PetscValidCharPointer(type,2); 142ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 143ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1441c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 14584df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 146ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1474cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 14884df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 14968ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15084df9cb4SJed Brown PetscFunctionReturn(0); 15184df9cb4SJed Brown } 15284df9cb4SJed Brown 153d0288e4fSLisandro Dalcin /*@C 154d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 155d0288e4fSLisandro Dalcin 156d0288e4fSLisandro Dalcin Not Collective 157d0288e4fSLisandro Dalcin 158d0288e4fSLisandro Dalcin Input Parameter: 159d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 160d0288e4fSLisandro Dalcin 161d0288e4fSLisandro Dalcin Output Parameter: 162d0288e4fSLisandro Dalcin . type - The name of TS adapter method 163d0288e4fSLisandro Dalcin 164d0288e4fSLisandro Dalcin Level: intermediate 165d0288e4fSLisandro Dalcin 166d0288e4fSLisandro Dalcin .seealso TSAdaptSetType() 167d0288e4fSLisandro Dalcin @*/ 168d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 169d0288e4fSLisandro Dalcin { 170d0288e4fSLisandro Dalcin PetscFunctionBegin; 171d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 172d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 173d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 174d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 175d0288e4fSLisandro Dalcin } 176d0288e4fSLisandro Dalcin 17784df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 17884df9cb4SJed Brown { 17984df9cb4SJed Brown PetscErrorCode ierr; 18084df9cb4SJed Brown 18184df9cb4SJed Brown PetscFunctionBegin; 1824782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 18384df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 18484df9cb4SJed Brown PetscFunctionReturn(0); 18584df9cb4SJed Brown } 18684df9cb4SJed Brown 187ad6bc421SBarry Smith /*@C 188ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 189ad6bc421SBarry Smith 190ad6bc421SBarry Smith Collective on PetscViewer 191ad6bc421SBarry Smith 192ad6bc421SBarry Smith Input Parameters: 193ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 194ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 195ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 196ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 197ad6bc421SBarry Smith 198ad6bc421SBarry Smith Level: intermediate 199ad6bc421SBarry Smith 200ad6bc421SBarry Smith Notes: 201ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 202ad6bc421SBarry Smith 203ad6bc421SBarry Smith Notes for advanced users: 204ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 205ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 206ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 207ad6bc421SBarry Smith format is 208ad6bc421SBarry Smith .vb 209ad6bc421SBarry Smith has not yet been determined 210ad6bc421SBarry Smith .ve 211ad6bc421SBarry Smith 212ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 213ad6bc421SBarry Smith @*/ 2144782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 215ad6bc421SBarry Smith { 216ad6bc421SBarry Smith PetscErrorCode ierr; 217ad6bc421SBarry Smith PetscBool isbinary; 218ad6bc421SBarry Smith char type[256]; 219ad6bc421SBarry Smith 220ad6bc421SBarry Smith PetscFunctionBegin; 2214782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 222ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 223ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 224ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 225ad6bc421SBarry Smith 226060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2274782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2284782b174SLisandro Dalcin if (adapt->ops->load) { 2294782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 230ad6bc421SBarry Smith } 231ad6bc421SBarry Smith PetscFunctionReturn(0); 232ad6bc421SBarry Smith } 233ad6bc421SBarry Smith 23484df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 23584df9cb4SJed Brown { 23684df9cb4SJed Brown PetscErrorCode ierr; 2371c167fc2SEmil Constantinescu PetscBool iascii,isbinary,isnone,isglee; 23884df9cb4SJed Brown 23984df9cb4SJed Brown PetscFunctionBegin; 2404782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2414782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2424782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2434782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 244251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 245ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 24684df9cb4SJed Brown if (iascii) { 247dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 2481917a363SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr); 2491c167fc2SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTGLEE,&isglee);CHKERRQ(ierr); 2501917a363SLisandro Dalcin if (!isnone) { 251bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 2521917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 2531917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 2541917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 2551917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 256bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 257bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 2581c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," maximum solution absolute value to be ignored %g\n",(double)adapt->ignore_max);CHKERRQ(ierr); 2591c167fc2SEmil Constantinescu } 2601c167fc2SEmil Constantinescu if (isglee) { 2611c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2621c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE uses local error control\n");CHKERRQ(ierr); 2631c167fc2SEmil Constantinescu } else { 2641c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE uses global error control\n");CHKERRQ(ierr); 2651c167fc2SEmil Constantinescu } 2661917a363SLisandro Dalcin } 26784df9cb4SJed Brown if (adapt->ops->view) { 26884df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 26984df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 27084df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27184df9cb4SJed Brown } 272ad6bc421SBarry Smith } else if (isbinary) { 273ad6bc421SBarry Smith char type[256]; 274ad6bc421SBarry Smith 275ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 276ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 277f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 278bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 279f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 280f2c2a1b9SBarry Smith } 28184df9cb4SJed Brown PetscFunctionReturn(0); 28284df9cb4SJed Brown } 28384df9cb4SJed Brown 284099af0a3SLisandro Dalcin /*@ 285099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 286099af0a3SLisandro Dalcin 287099af0a3SLisandro Dalcin Collective on TS 288099af0a3SLisandro Dalcin 289099af0a3SLisandro Dalcin Input Parameter: 290099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 291099af0a3SLisandro Dalcin 292099af0a3SLisandro Dalcin Level: developer 293099af0a3SLisandro Dalcin 294099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 295099af0a3SLisandro Dalcin @*/ 296099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 297099af0a3SLisandro Dalcin { 298099af0a3SLisandro Dalcin PetscErrorCode ierr; 299099af0a3SLisandro Dalcin 300099af0a3SLisandro Dalcin PetscFunctionBegin; 301099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 302099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 303099af0a3SLisandro Dalcin PetscFunctionReturn(0); 304099af0a3SLisandro Dalcin } 305099af0a3SLisandro Dalcin 30684df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 30784df9cb4SJed Brown { 30884df9cb4SJed Brown PetscErrorCode ierr; 30984df9cb4SJed Brown 31084df9cb4SJed Brown PetscFunctionBegin; 31184df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 31284df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3134782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 314099af0a3SLisandro Dalcin 315099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 316099af0a3SLisandro Dalcin 31784df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3181c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 31984df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 32084df9cb4SJed Brown PetscFunctionReturn(0); 32184df9cb4SJed Brown } 32284df9cb4SJed Brown 3231c3436cfSJed Brown /*@ 3241c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3251c3436cfSJed Brown 3261c3436cfSJed Brown Collective on TSAdapt 3271c3436cfSJed Brown 3281c3436cfSJed Brown Input Arguments: 3291c3436cfSJed Brown + adapt - adaptive controller context 3301c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3311c3436cfSJed Brown 332bf997491SLisandro Dalcin Options Database Keys: 333ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 334bf997491SLisandro Dalcin 3351c3436cfSJed Brown Level: intermediate 3361c3436cfSJed Brown 3371c3436cfSJed Brown .seealso: TSAdaptChoose() 3381c3436cfSJed Brown @*/ 3391c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3401c3436cfSJed Brown { 3411c3436cfSJed Brown PetscErrorCode ierr; 3421c3436cfSJed Brown 3431c3436cfSJed Brown PetscFunctionBegin; 3444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3454782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3461c3436cfSJed Brown if (flg) { 347ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3481c3436cfSJed Brown } else { 3491c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3501c3436cfSJed Brown } 3511c3436cfSJed Brown PetscFunctionReturn(0); 3521c3436cfSJed Brown } 3531c3436cfSJed Brown 3540873808bSJed Brown /*@C 355bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3560873808bSJed Brown 3570873808bSJed Brown Logically collective on TSAdapt 3580873808bSJed Brown 3590873808bSJed Brown Input Arguments: 3600873808bSJed Brown + adapt - adaptive controller context 3610873808bSJed Brown - func - stage check function 3620873808bSJed Brown 3630873808bSJed Brown Arguments of func: 3640873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3650873808bSJed Brown 3660873808bSJed Brown + adapt - adaptive controller context 3670873808bSJed Brown . ts - time stepping context 3680873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3690873808bSJed Brown 3700873808bSJed Brown Level: advanced 3710873808bSJed Brown 3720873808bSJed Brown .seealso: TSAdaptChoose() 3730873808bSJed Brown @*/ 374b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3750873808bSJed Brown { 3760873808bSJed Brown PetscFunctionBegin; 3770873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 37868ae941aSLisandro Dalcin adapt->checkstage = func; 3790873808bSJed Brown PetscFunctionReturn(0); 3800873808bSJed Brown } 3810873808bSJed Brown 3821c3436cfSJed Brown /*@ 383bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 384bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 385bf997491SLisandro Dalcin 3861917a363SLisandro Dalcin Logically collective on TSAdapt 387bf997491SLisandro Dalcin 388bf997491SLisandro Dalcin Input Arguments: 389bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 390bf997491SLisandro Dalcin - flag - whether to always accept steps 391bf997491SLisandro Dalcin 392bf997491SLisandro Dalcin Options Database Keys: 393ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 394bf997491SLisandro Dalcin 395bf997491SLisandro Dalcin Level: intermediate 396bf997491SLisandro Dalcin 397bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 398bf997491SLisandro Dalcin @*/ 399bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 400bf997491SLisandro Dalcin { 401bf997491SLisandro Dalcin PetscFunctionBegin; 402bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 403bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 404bf997491SLisandro Dalcin adapt->always_accept = flag; 405bf997491SLisandro Dalcin PetscFunctionReturn(0); 406bf997491SLisandro Dalcin } 407bf997491SLisandro Dalcin 408bf997491SLisandro Dalcin /*@ 4091917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4101c3436cfSJed Brown 4111917a363SLisandro Dalcin Logically collective on TSAdapt 4121917a363SLisandro Dalcin 4131917a363SLisandro Dalcin Input Arguments: 4141917a363SLisandro Dalcin + adapt - adaptive controller context 4151917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 416ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 4171917a363SLisandro Dalcin 4181917a363SLisandro Dalcin Options Database Keys: 419ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 420ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4211917a363SLisandro Dalcin 4221917a363SLisandro Dalcin Level: intermediate 4231917a363SLisandro Dalcin 4241917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4251917a363SLisandro Dalcin @*/ 4261917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4271917a363SLisandro Dalcin { 4281917a363SLisandro Dalcin PetscFunctionBegin; 4291917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4301917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4311917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4321917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4331917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety); 4341917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT && reject_safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be non negative",(double)reject_safety); 4351917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT && reject_safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be less than one",(double)reject_safety); 4361917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4371917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4381917a363SLisandro Dalcin PetscFunctionReturn(0); 4391917a363SLisandro Dalcin } 4401917a363SLisandro Dalcin 4411917a363SLisandro Dalcin /*@ 4421917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4431917a363SLisandro Dalcin 4441917a363SLisandro Dalcin Not Collective 4451917a363SLisandro Dalcin 4461917a363SLisandro Dalcin Input Arguments: 4471917a363SLisandro Dalcin . adapt - adaptive controller context 4481917a363SLisandro Dalcin 44901d2d390SJose E. Roman Output Arguments: 4501917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4511917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4521917a363SLisandro Dalcin 4531917a363SLisandro Dalcin Level: intermediate 4541917a363SLisandro Dalcin 4551917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4561917a363SLisandro Dalcin @*/ 4571917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4581917a363SLisandro Dalcin { 4591917a363SLisandro Dalcin PetscFunctionBegin; 4601917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4611917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4621917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4631917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4641917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4651917a363SLisandro Dalcin PetscFunctionReturn(0); 4661917a363SLisandro Dalcin } 4671917a363SLisandro Dalcin 4681917a363SLisandro Dalcin /*@ 46976cddca1SEmil 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. 47076cddca1SEmil Constantinescu 47176cddca1SEmil Constantinescu Logically collective on TSAdapt 47276cddca1SEmil Constantinescu 47376cddca1SEmil Constantinescu Input Arguments: 47476cddca1SEmil Constantinescu + adapt - adaptive controller context 47576cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 47676cddca1SEmil Constantinescu 47776cddca1SEmil Constantinescu Options Database Keys: 47876cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 47976cddca1SEmil Constantinescu 48076cddca1SEmil Constantinescu Level: intermediate 48176cddca1SEmil Constantinescu 48276cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptGetMaxIgnore(), TSAdaptChoose() 48376cddca1SEmil Constantinescu @*/ 48476cddca1SEmil Constantinescu PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore) 48576cddca1SEmil Constantinescu { 48676cddca1SEmil Constantinescu PetscFunctionBegin; 48776cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 48876cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt,max_ignore,2); 48976cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 49076cddca1SEmil Constantinescu PetscFunctionReturn(0); 49176cddca1SEmil Constantinescu } 49276cddca1SEmil Constantinescu 49376cddca1SEmil Constantinescu /*@ 49476cddca1SEmil 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). 49576cddca1SEmil Constantinescu 49676cddca1SEmil Constantinescu Not Collective 49776cddca1SEmil Constantinescu 49876cddca1SEmil Constantinescu Input Arguments: 49976cddca1SEmil Constantinescu . adapt - adaptive controller context 50076cddca1SEmil Constantinescu 50101d2d390SJose E. Roman Output Arguments: 50276cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 50376cddca1SEmil Constantinescu 50476cddca1SEmil Constantinescu Level: intermediate 50576cddca1SEmil Constantinescu 50676cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptSetMaxIgnore(), TSAdaptChoose() 50776cddca1SEmil Constantinescu @*/ 50876cddca1SEmil Constantinescu PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore) 50976cddca1SEmil Constantinescu { 51076cddca1SEmil Constantinescu PetscFunctionBegin; 51176cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 51276cddca1SEmil Constantinescu PetscValidRealPointer(max_ignore,2); 51376cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 51476cddca1SEmil Constantinescu PetscFunctionReturn(0); 51576cddca1SEmil Constantinescu } 51676cddca1SEmil Constantinescu 51776cddca1SEmil Constantinescu /*@ 5181917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 5191917a363SLisandro Dalcin 5201917a363SLisandro Dalcin Logically collective on TSAdapt 5211917a363SLisandro Dalcin 5221917a363SLisandro Dalcin Input Arguments: 5231917a363SLisandro Dalcin + adapt - adaptive controller context 5241917a363SLisandro Dalcin . low - admissible decrease factor 525ec18b7bcSLisandro Dalcin - high - admissible increase factor 5261917a363SLisandro Dalcin 5271917a363SLisandro Dalcin Options Database Keys: 528ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5291917a363SLisandro Dalcin 5301917a363SLisandro Dalcin Level: intermediate 5311917a363SLisandro Dalcin 53262c23b28SMark .seealso: TSAdaptChoose(), TSAdaptGetClip(), TSAdaptSetScaleSolveFailed() 5331917a363SLisandro Dalcin @*/ 5341917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 5351917a363SLisandro Dalcin { 5361917a363SLisandro Dalcin PetscFunctionBegin; 5371917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5381917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 5391917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 5401917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 5411917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 5422479783cSJose E. Roman if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be greater than one",(double)high); 5431917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 5441917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 5451917a363SLisandro Dalcin PetscFunctionReturn(0); 5461917a363SLisandro Dalcin } 5471917a363SLisandro Dalcin 5481917a363SLisandro Dalcin /*@ 5491917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5501917a363SLisandro Dalcin 5511917a363SLisandro Dalcin Not Collective 5521917a363SLisandro Dalcin 5531917a363SLisandro Dalcin Input Arguments: 5541917a363SLisandro Dalcin . adapt - adaptive controller context 5551917a363SLisandro Dalcin 55601d2d390SJose E. Roman Output Arguments: 5571917a363SLisandro Dalcin + low - optional, admissible decrease factor 5581917a363SLisandro Dalcin - high - optional, admissible increase factor 5591917a363SLisandro Dalcin 5601917a363SLisandro Dalcin Level: intermediate 5611917a363SLisandro Dalcin 56262c23b28SMark .seealso: TSAdaptChoose(), TSAdaptSetClip(), TSAdaptSetScaleSolveFailed() 5631917a363SLisandro Dalcin @*/ 5641917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5651917a363SLisandro Dalcin { 5661917a363SLisandro Dalcin PetscFunctionBegin; 5671917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5681917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5691917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5701917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5711917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5721917a363SLisandro Dalcin PetscFunctionReturn(0); 5731917a363SLisandro Dalcin } 5741917a363SLisandro Dalcin 5751917a363SLisandro Dalcin /*@ 57662c23b28SMark TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 57762c23b28SMark 57862c23b28SMark Logically collective on TSAdapt 57962c23b28SMark 58062c23b28SMark Input Arguments: 58162c23b28SMark + adapt - adaptive controller context 582ee300463SSatish Balay - scale - scale 58362c23b28SMark 58462c23b28SMark Options Database Keys: 58562c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 58662c23b28SMark 58762c23b28SMark Level: intermediate 58862c23b28SMark 58962c23b28SMark .seealso: TSAdaptChoose(), TSAdaptGetScaleSolveFailed(), TSAdaptGetClip() 59062c23b28SMark @*/ 59162c23b28SMark PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt,PetscReal scale) 59262c23b28SMark { 59362c23b28SMark PetscFunctionBegin; 59462c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 59562c23b28SMark PetscValidLogicalCollectiveReal(adapt,scale,2); 59662c23b28SMark if (scale != PETSC_DEFAULT && scale <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be positive",(double)scale); 59762c23b28SMark if (scale != PETSC_DEFAULT && scale > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be less than one",(double)scale); 59862c23b28SMark if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 59962c23b28SMark PetscFunctionReturn(0); 60062c23b28SMark } 60162c23b28SMark 60262c23b28SMark /*@ 60362c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 60462c23b28SMark 60562c23b28SMark Not Collective 60662c23b28SMark 60762c23b28SMark Input Arguments: 60862c23b28SMark . adapt - adaptive controller context 60962c23b28SMark 61001d2d390SJose E. Roman Output Arguments: 611ee300463SSatish Balay . scale - scale factor 61262c23b28SMark 61362c23b28SMark Level: intermediate 61462c23b28SMark 61562c23b28SMark .seealso: TSAdaptChoose(), TSAdaptSetScaleSolveFailed(), TSAdaptSetClip() 61662c23b28SMark @*/ 61762c23b28SMark PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt,PetscReal *scale) 61862c23b28SMark { 61962c23b28SMark PetscFunctionBegin; 62062c23b28SMark PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 62162c23b28SMark if (scale) PetscValidRealPointer(scale,2); 62262c23b28SMark if (scale) *scale = adapt->scale_solve_failed; 62362c23b28SMark PetscFunctionReturn(0); 62462c23b28SMark } 62562c23b28SMark 62662c23b28SMark /*@ 6271917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 6281917a363SLisandro Dalcin 6291917a363SLisandro Dalcin Logically collective on TSAdapt 6301c3436cfSJed Brown 6311c3436cfSJed Brown Input Arguments: 632552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6331c3436cfSJed Brown . hmin - minimum time step 6341c3436cfSJed Brown - hmax - maximum time step 6351c3436cfSJed Brown 6361c3436cfSJed Brown Options Database Keys: 637ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 638ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 6391c3436cfSJed Brown 6401c3436cfSJed Brown Level: intermediate 6411c3436cfSJed Brown 6421917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 6431c3436cfSJed Brown @*/ 6441c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 6451c3436cfSJed Brown { 6461c3436cfSJed Brown 6471c3436cfSJed Brown PetscFunctionBegin; 6484782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 649b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 650b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 651b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 652b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 653b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 654b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 655b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 656b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 6572479783cSJose E. Roman if (hmax <= hmin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Maximum time step %g must greater than minimum time step %g",(double)hmax,(double)hmin); 6581917a363SLisandro Dalcin PetscFunctionReturn(0); 6591917a363SLisandro Dalcin } 6601917a363SLisandro Dalcin 6611917a363SLisandro Dalcin /*@ 6621917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 6631917a363SLisandro Dalcin 6641917a363SLisandro Dalcin Not Collective 6651917a363SLisandro Dalcin 6661917a363SLisandro Dalcin Input Arguments: 6671917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6681917a363SLisandro Dalcin 6691917a363SLisandro Dalcin Output Arguments: 6701917a363SLisandro Dalcin + hmin - minimum time step 6711917a363SLisandro Dalcin - hmax - maximum time step 6721917a363SLisandro Dalcin 6731917a363SLisandro Dalcin Level: intermediate 6741917a363SLisandro Dalcin 6751917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 6761917a363SLisandro Dalcin @*/ 6771917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 6781917a363SLisandro Dalcin { 6791917a363SLisandro Dalcin 6801917a363SLisandro Dalcin PetscFunctionBegin; 6811917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6821917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 6831917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 6841917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6851917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6861c3436cfSJed Brown PetscFunctionReturn(0); 6871c3436cfSJed Brown } 6881c3436cfSJed Brown 689e55864a3SBarry Smith /* 69084df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 69184df9cb4SJed Brown 69284df9cb4SJed Brown Collective on TSAdapt 69384df9cb4SJed Brown 69484df9cb4SJed Brown Input Parameter: 69584df9cb4SJed Brown . adapt - the TSAdapt context 69684df9cb4SJed Brown 69784df9cb4SJed Brown Options Database Keys: 6981917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 699bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 7001917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 7011917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 7021917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 70323de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 70423de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 70523de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 706de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 707de50f1caSBarry 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 70884df9cb4SJed Brown 70984df9cb4SJed Brown Level: advanced 71084df9cb4SJed Brown 71184df9cb4SJed Brown Notes: 71284df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 71384df9cb4SJed Brown 7141917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 71562c23b28SMark TSAdaptSetClip(), TSAdaptSetScaleSolveFailed(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 716e55864a3SBarry Smith */ 7174416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 71884df9cb4SJed Brown { 71984df9cb4SJed Brown PetscErrorCode ierr; 72084df9cb4SJed Brown char type[256] = TSADAPTBASIC; 72162c23b28SMark PetscReal safety,reject_safety,clip[2],scale,hmin,hmax; 7221c3436cfSJed Brown PetscBool set,flg; 7231917a363SLisandro Dalcin PetscInt two; 72484df9cb4SJed Brown 72584df9cb4SJed Brown PetscFunctionBegin; 726064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,2); 72784df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 7281566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 729e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 73023de1d84SBarry Smith ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 73184df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 73284df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 73384df9cb4SJed Brown } 7341917a363SLisandro Dalcin 735bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 736bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 7371917a363SLisandro Dalcin 7381917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 7391917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 7401917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg);CHKERRQ(ierr); 7411917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 7421917a363SLisandro Dalcin 7431917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 7441917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 7451917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 7461917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 7471917a363SLisandro Dalcin 7481917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 7491917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 7501917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 751bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 7521917a363SLisandro Dalcin 753d580f011SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set);CHKERRQ(ierr); 754d580f011SEmil Constantinescu ierr = PetscOptionsBool("-ts_adapt_glee_use_local","GLEE adaptor uses local error estimation for step control","",adapt->glee_use_local,&adapt->glee_use_local,&set);CHKERRQ(ierr); 755d580f011SEmil Constantinescu 75662c23b28SMark ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","TSAdaptSetScaleSolveFailed",adapt->scale_solve_failed,&scale,&set);CHKERRQ(ierr); 75762c23b28SMark if (set) {ierr = TSAdaptSetScaleSolveFailed(adapt,scale);CHKERRQ(ierr);} 7581917a363SLisandro Dalcin 7597619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 7607619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 7611917a363SLisandro Dalcin 762de50f1caSBarry Smith ierr = 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);CHKERRQ(ierr); 763de50f1caSBarry Smith 7641917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 7651917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 7661917a363SLisandro Dalcin 767e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 76884df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 76984df9cb4SJed Brown PetscFunctionReturn(0); 77084df9cb4SJed Brown } 77184df9cb4SJed Brown 77284df9cb4SJed Brown /*@ 77384df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 77484df9cb4SJed Brown 7751917a363SLisandro Dalcin Logically collective on TSAdapt 77684df9cb4SJed Brown 77784df9cb4SJed Brown Input Argument: 77884df9cb4SJed Brown . adapt - adaptive controller 77984df9cb4SJed Brown 78084df9cb4SJed Brown Level: developer 78184df9cb4SJed Brown 78284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 78384df9cb4SJed Brown @*/ 78484df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 78584df9cb4SJed Brown { 78684df9cb4SJed Brown PetscErrorCode ierr; 78784df9cb4SJed Brown 78884df9cb4SJed Brown PetscFunctionBegin; 7894782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 79084df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 79184df9cb4SJed Brown PetscFunctionReturn(0); 79284df9cb4SJed Brown } 79384df9cb4SJed Brown 79484df9cb4SJed Brown /*@C 79584df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 79684df9cb4SJed Brown 7971917a363SLisandro Dalcin Logically collective on TSAdapt 79884df9cb4SJed Brown 79984df9cb4SJed Brown Input Arguments: 800552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 80184df9cb4SJed Brown . name - name of the candidate scheme to add 80284df9cb4SJed Brown . order - order of the candidate scheme 80384df9cb4SJed Brown . stageorder - stage order of the candidate scheme 8048d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 80584df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 80684df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 80784df9cb4SJed Brown 80884df9cb4SJed Brown Note: 80984df9cb4SJed Brown This routine is not available in Fortran. 81084df9cb4SJed Brown 81184df9cb4SJed Brown Level: developer 81284df9cb4SJed Brown 81384df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 81484df9cb4SJed Brown @*/ 8158d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 81684df9cb4SJed Brown { 81784df9cb4SJed Brown PetscInt c; 81884df9cb4SJed Brown 81984df9cb4SJed Brown PetscFunctionBegin; 82084df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 821ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 82284df9cb4SJed Brown if (inuse) { 823ce94432eSBarry Smith if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 82484df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 82584df9cb4SJed Brown } 8261c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 8271c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 828bbd56ea5SKarl Rupp 82984df9cb4SJed Brown adapt->candidates.name[c] = name; 83084df9cb4SJed Brown adapt->candidates.order[c] = order; 83184df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 8328d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 83384df9cb4SJed Brown adapt->candidates.cost[c] = cost; 83484df9cb4SJed Brown adapt->candidates.n++; 83584df9cb4SJed Brown PetscFunctionReturn(0); 83684df9cb4SJed Brown } 83784df9cb4SJed Brown 8388d59e960SJed Brown /*@C 8398d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8408d59e960SJed Brown 8418d59e960SJed Brown Not Collective 8428d59e960SJed Brown 8438d59e960SJed Brown Input Arguments: 8448d59e960SJed Brown . adapt - time step adaptivity context 8458d59e960SJed Brown 8468d59e960SJed Brown Output Arguments: 8478d59e960SJed Brown + n - number of candidate schemes, always at least 1 8488d59e960SJed Brown . order - the order of each candidate scheme 8498d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8508d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8518d59e960SJed Brown - cost - the relative cost of each scheme 8528d59e960SJed Brown 8538d59e960SJed Brown Level: developer 8548d59e960SJed Brown 8558d59e960SJed Brown Note: 8568d59e960SJed Brown The current scheme is always returned in the first slot 8578d59e960SJed Brown 8588d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 8598d59e960SJed Brown @*/ 8608d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 8618d59e960SJed Brown { 8628d59e960SJed Brown PetscFunctionBegin; 8638d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8648d59e960SJed Brown if (n) *n = adapt->candidates.n; 8658d59e960SJed Brown if (order) *order = adapt->candidates.order; 8668d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8678d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8688d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8698d59e960SJed Brown PetscFunctionReturn(0); 8708d59e960SJed Brown } 8718d59e960SJed Brown 87284df9cb4SJed Brown /*@C 87384df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 87484df9cb4SJed Brown 8751917a363SLisandro Dalcin Collective on TSAdapt 87684df9cb4SJed Brown 87784df9cb4SJed Brown Input Arguments: 87884df9cb4SJed Brown + adapt - adaptive contoller 87984df9cb4SJed Brown - h - current step size 88084df9cb4SJed Brown 88184df9cb4SJed Brown Output Arguments: 8821566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 88384df9cb4SJed Brown . next_h - step size to use for the next step 88484df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 88584df9cb4SJed Brown 8861c3436cfSJed Brown Note: 8871c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 8881c3436cfSJed Brown being retried after an initial rejection. 8891c3436cfSJed Brown 89084df9cb4SJed Brown Level: developer 89184df9cb4SJed Brown 89284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 89384df9cb4SJed Brown @*/ 89484df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 89584df9cb4SJed Brown { 89684df9cb4SJed Brown PetscErrorCode ierr; 8971566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8981566a47fSLisandro Dalcin PetscInt scheme = 0; 8990b99f514SJed Brown PetscReal wlte = -1.0; 9005e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 9015e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 90284df9cb4SJed Brown 90384df9cb4SJed Brown PetscFunctionBegin; 90484df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 90584df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 9061566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 90784df9cb4SJed Brown PetscValidPointer(next_h,5); 908064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,6); 9091566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 9101566a47fSLisandro Dalcin 9111566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 9121566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 9131566a47fSLisandro Dalcin *next_h = h; 9141566a47fSLisandro Dalcin *accept = PETSC_TRUE; 9151566a47fSLisandro Dalcin PetscFunctionReturn(0); 9161566a47fSLisandro Dalcin } 9171566a47fSLisandro Dalcin 9185e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 9191566a47fSLisandro Dalcin if (scheme < 0 || (ncandidates > 0 && scheme >= ncandidates)) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",scheme,ncandidates-1); 9201566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 9211566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 9221566a47fSLisandro Dalcin 92349354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 92436b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 92536b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 92636b54a69SLisandro Dalcin PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 927fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 928fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 92936b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 930fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 93136b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 93249354f04SShri Abhyankar } 9331c3436cfSJed Brown 9341c3436cfSJed Brown if (adapt->monitor) { 9351566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 9361c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9370b99f514SJed Brown if (wlte < 0) { 9380dea241dSBarry Smith ierr = 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);CHKERRQ(ierr); 9390b99f514SJed Brown } else { 9400dea241dSBarry Smith ierr = 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);CHKERRQ(ierr); 9410b99f514SJed Brown } 9421c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9431c3436cfSJed Brown } 94484df9cb4SJed Brown PetscFunctionReturn(0); 94584df9cb4SJed Brown } 94684df9cb4SJed Brown 94797335746SJed Brown /*@ 948de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 949de50f1caSBarry Smith before increasing the time step. 950de50f1caSBarry Smith 951de50f1caSBarry Smith Logicially Collective on TSAdapt 952de50f1caSBarry Smith 953de50f1caSBarry Smith Input Arguments: 954de50f1caSBarry Smith + adapt - adaptive controller context 955de50f1caSBarry Smith - cnt - the number of timesteps 956de50f1caSBarry Smith 957de50f1caSBarry Smith Options Database Key: 958de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 959de50f1caSBarry Smith 960de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 961de50f1caSBarry Smith The successful use of this option is problem dependent 962de50f1caSBarry Smith 963de50f1caSBarry Smith Developer Note: there is no theory to support this option 964de50f1caSBarry Smith 965de50f1caSBarry Smith Level: advanced 966de50f1caSBarry Smith 967de50f1caSBarry Smith .seealso: 968de50f1caSBarry Smith @*/ 969de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 970de50f1caSBarry Smith { 971de50f1caSBarry Smith PetscFunctionBegin; 972de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 973de50f1caSBarry Smith PetscFunctionReturn(0); 974de50f1caSBarry Smith } 975de50f1caSBarry Smith 976de50f1caSBarry Smith /*@ 9776bc98fa9SBarry 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) 97897335746SJed Brown 9791917a363SLisandro Dalcin Collective on TSAdapt 98097335746SJed Brown 98197335746SJed Brown Input Arguments: 98297335746SJed Brown + adapt - adaptive controller context 983b295832fSPierre Barbier de Reuille . ts - time stepper 984b295832fSPierre Barbier de Reuille . t - Current simulation time 985b295832fSPierre Barbier de Reuille - Y - Current solution vector 98697335746SJed Brown 98797335746SJed Brown Output Arguments: 98897335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 98997335746SJed Brown 99097335746SJed Brown Level: developer 99197335746SJed Brown 99297335746SJed Brown .seealso: 99397335746SJed Brown @*/ 994b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 99597335746SJed Brown { 99697335746SJed Brown PetscErrorCode ierr; 9971566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 99897335746SJed Brown 99997335746SJed Brown PetscFunctionBegin; 10004782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 10014782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1002064a246eSJacob Faibussowitsch PetscValidBoolPointer(accept,5); 10031566a47fSLisandro Dalcin 10041566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 100597335746SJed Brown if (snesreason < 0) { 100697335746SJed Brown *accept = PETSC_FALSE; 10076de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 100897335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 100948c19aefSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 101097335746SJed Brown if (adapt->monitor) { 101197335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 10120dea241dSBarry Smith ierr = 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);CHKERRQ(ierr); 101397335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 101497335746SJed Brown } 1015cb9d8021SPierre Barbier de Reuille } 1016cb9d8021SPierre Barbier de Reuille } else { 10171566a47fSLisandro Dalcin *accept = PETSC_TRUE; 1018b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 1019cb9d8021SPierre Barbier de Reuille if (*accept && adapt->checkstage) { 1020b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 10216bc98fa9SBarry Smith if (!*accept) { 10226bc98fa9SBarry Smith ierr = PetscInfo1(ts,"Step=%D, solution rejected by user function provided by TSSetFunctionDomainError()\n",ts->steps);CHKERRQ(ierr); 10236bc98fa9SBarry Smith if (adapt->monitor) { 10246bc98fa9SBarry Smith ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 10256bc98fa9SBarry Smith ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected by user function provided by TSSetFunctionDomainError()\n",((PetscObject)adapt)->type_name,ts->steps);CHKERRQ(ierr); 10266bc98fa9SBarry Smith ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 10276bc98fa9SBarry Smith } 10286bc98fa9SBarry Smith } 1029cb9d8021SPierre Barbier de Reuille } 1030cb9d8021SPierre Barbier de Reuille } 1031cb9d8021SPierre Barbier de Reuille 10321566a47fSLisandro Dalcin if (!(*accept) && !ts->reason) { 10331566a47fSLisandro Dalcin PetscReal dt,new_dt; 10341566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 1035cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 103697335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 1037de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 103897335746SJed Brown if (adapt->monitor) { 103997335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 1040e6d0a238SBarry Smith ierr = 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);CHKERRQ(ierr); 104197335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 104297335746SJed Brown } 104397335746SJed Brown } 104497335746SJed Brown PetscFunctionReturn(0); 104597335746SJed Brown } 104697335746SJed Brown 104784df9cb4SJed Brown /*@ 104884df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 104984df9cb4SJed Brown 1050d083f849SBarry Smith Collective 105184df9cb4SJed Brown 105284df9cb4SJed Brown Input Parameter: 105384df9cb4SJed Brown . comm - The communicator 105484df9cb4SJed Brown 105584df9cb4SJed Brown Output Parameter: 105684df9cb4SJed Brown . adapt - new TSAdapt object 105784df9cb4SJed Brown 105884df9cb4SJed Brown Level: developer 105984df9cb4SJed Brown 106084df9cb4SJed Brown Notes: 106184df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 106284df9cb4SJed Brown 1063552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 106484df9cb4SJed Brown @*/ 106584df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 106684df9cb4SJed Brown { 106784df9cb4SJed Brown PetscErrorCode ierr; 106884df9cb4SJed Brown TSAdapt adapt; 106984df9cb4SJed Brown 107084df9cb4SJed Brown PetscFunctionBegin; 1071064a246eSJacob Faibussowitsch PetscValidPointer(inadapt,2); 10723b3bcf4cSLisandro Dalcin *inadapt = NULL; 10733b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 10743b3bcf4cSLisandro Dalcin 107573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 10761c3436cfSJed Brown 1077bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 10781917a363SLisandro Dalcin adapt->safety = 0.9; 10791917a363SLisandro Dalcin adapt->reject_safety = 0.5; 10801917a363SLisandro Dalcin adapt->clip[0] = 0.1; 10811917a363SLisandro Dalcin adapt->clip[1] = 10.; 10821c3436cfSJed Brown adapt->dt_min = 1e-20; 10831917a363SLisandro Dalcin adapt->dt_max = 1e+20; 10841c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1085d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 108697335746SJed Brown adapt->scale_solve_failed = 0.25; 1087fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1088fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1089fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1090fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 10917619abb3SShri adapt->wnormtype = NORM_2; 1092de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 10931917a363SLisandro Dalcin 109484df9cb4SJed Brown *inadapt = adapt; 109584df9cb4SJed Brown PetscFunctionReturn(0); 109684df9cb4SJed Brown } 1097