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 .keywords: TSAdapt, register 421c84c290SBarry Smith 431c84c290SBarry Smith .seealso: TSAdaptRegisterAll() 4484df9cb4SJed Brown @*/ 45bdf89e91SBarry Smith PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 4684df9cb4SJed Brown { 4784df9cb4SJed Brown PetscErrorCode ierr; 4884df9cb4SJed Brown 4984df9cb4SJed Brown PetscFunctionBegin; 501d36bdfdSBarry Smith ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 51a240a19fSJed Brown ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 5284df9cb4SJed Brown PetscFunctionReturn(0); 5384df9cb4SJed Brown } 5484df9cb4SJed Brown 5584df9cb4SJed Brown /*@C 5684df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5784df9cb4SJed Brown 5884df9cb4SJed Brown Not Collective 5984df9cb4SJed Brown 6084df9cb4SJed Brown Level: advanced 6184df9cb4SJed Brown 6284df9cb4SJed Brown .keywords: TSAdapt, register, all 6384df9cb4SJed Brown 6484df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6584df9cb4SJed Brown @*/ 66607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6784df9cb4SJed Brown { 6884df9cb4SJed Brown PetscErrorCode ierr; 6984df9cb4SJed Brown 7084df9cb4SJed Brown PetscFunctionBegin; 710f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 720f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 73bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 741566a47fSLisandro Dalcin ierr = TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic);CHKERRQ(ierr); 75cb7ab186SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP);CHKERRQ(ierr); 76bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 771917a363SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE);CHKERRQ(ierr); 78d949e4a4SStefano Zampini ierr = TSAdaptRegister(TSADAPTHISTORY,TSAdaptCreate_History);CHKERRQ(ierr); 7984df9cb4SJed Brown PetscFunctionReturn(0); 8084df9cb4SJed Brown } 8184df9cb4SJed Brown 8284df9cb4SJed Brown /*@C 83560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 8484df9cb4SJed Brown called from PetscFinalize(). 8584df9cb4SJed Brown 8684df9cb4SJed Brown Level: developer 8784df9cb4SJed Brown 8884df9cb4SJed Brown .keywords: Petsc, destroy, package 8984df9cb4SJed Brown .seealso: PetscFinalize() 9084df9cb4SJed Brown @*/ 9184df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 9284df9cb4SJed Brown { 9337e93019SBarry Smith PetscErrorCode ierr; 9437e93019SBarry Smith 9584df9cb4SJed Brown PetscFunctionBegin; 9637e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9784df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9884df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9984df9cb4SJed Brown PetscFunctionReturn(0); 10084df9cb4SJed Brown } 10184df9cb4SJed Brown 10284df9cb4SJed Brown /*@C 10384df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 1048a690491SBarry Smith called from TSInitializePackage(). 10584df9cb4SJed Brown 10684df9cb4SJed Brown Level: developer 10784df9cb4SJed Brown 10884df9cb4SJed Brown .keywords: TSAdapt, initialize, package 10984df9cb4SJed Brown .seealso: PetscInitialize() 11084df9cb4SJed Brown @*/ 111607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11284df9cb4SJed Brown { 11384df9cb4SJed Brown PetscErrorCode ierr; 11484df9cb4SJed Brown 11584df9cb4SJed Brown PetscFunctionBegin; 11684df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11784df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11884df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 119607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 12084df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 12184df9cb4SJed Brown PetscFunctionReturn(0); 12284df9cb4SJed Brown } 12384df9cb4SJed Brown 1247eef6055SBarry Smith /*@C 1257eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1267eef6055SBarry Smith 1277eef6055SBarry Smith Logicially Collective on TSAdapt 1287eef6055SBarry Smith 1297eef6055SBarry Smith Input Parameter: 130d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1317eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1327eef6055SBarry Smith 1337eef6055SBarry Smith Options Database: 134ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1357eef6055SBarry Smith 1367eef6055SBarry Smith Level: intermediate 1377eef6055SBarry Smith 1387eef6055SBarry Smith .keywords: TSAdapt, create 1397eef6055SBarry Smith 1407eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1417eef6055SBarry Smith @*/ 14219fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 14384df9cb4SJed Brown { 144ef749922SLisandro Dalcin PetscBool match; 14584df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 14684df9cb4SJed Brown 14784df9cb4SJed Brown PetscFunctionBegin; 1484782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 149b92453a8SLisandro Dalcin PetscValidCharPointer(type,2); 150ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 151ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1521c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 15384df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 154ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1554cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 15684df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15768ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15884df9cb4SJed Brown PetscFunctionReturn(0); 15984df9cb4SJed Brown } 16084df9cb4SJed Brown 161d0288e4fSLisandro Dalcin /*@C 162d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 163d0288e4fSLisandro Dalcin 164d0288e4fSLisandro Dalcin Not Collective 165d0288e4fSLisandro Dalcin 166d0288e4fSLisandro Dalcin Input Parameter: 167d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 168d0288e4fSLisandro Dalcin 169d0288e4fSLisandro Dalcin Output Parameter: 170d0288e4fSLisandro Dalcin . type - The name of TS adapter method 171d0288e4fSLisandro Dalcin 172d0288e4fSLisandro Dalcin Level: intermediate 173d0288e4fSLisandro Dalcin 174d0288e4fSLisandro Dalcin .keywords: TSAdapt, get, type 175d0288e4fSLisandro Dalcin .seealso TSAdaptSetType() 176d0288e4fSLisandro Dalcin @*/ 177d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 178d0288e4fSLisandro Dalcin { 179d0288e4fSLisandro Dalcin PetscFunctionBegin; 180d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 181d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 182d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 183d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 184d0288e4fSLisandro Dalcin } 185d0288e4fSLisandro Dalcin 18684df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 18784df9cb4SJed Brown { 18884df9cb4SJed Brown PetscErrorCode ierr; 18984df9cb4SJed Brown 19084df9cb4SJed Brown PetscFunctionBegin; 1914782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 19284df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 19384df9cb4SJed Brown PetscFunctionReturn(0); 19484df9cb4SJed Brown } 19584df9cb4SJed Brown 196ad6bc421SBarry Smith /*@C 197ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 198ad6bc421SBarry Smith 199ad6bc421SBarry Smith Collective on PetscViewer 200ad6bc421SBarry Smith 201ad6bc421SBarry Smith Input Parameters: 202ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 203ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 204ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 205ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 206ad6bc421SBarry Smith 207ad6bc421SBarry Smith Level: intermediate 208ad6bc421SBarry Smith 209ad6bc421SBarry Smith Notes: 210ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 211ad6bc421SBarry Smith 212ad6bc421SBarry Smith Notes for advanced users: 213ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 214ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 215ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 216ad6bc421SBarry Smith format is 217ad6bc421SBarry Smith .vb 218ad6bc421SBarry Smith has not yet been determined 219ad6bc421SBarry Smith .ve 220ad6bc421SBarry Smith 221ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 222ad6bc421SBarry Smith @*/ 2234782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 224ad6bc421SBarry Smith { 225ad6bc421SBarry Smith PetscErrorCode ierr; 226ad6bc421SBarry Smith PetscBool isbinary; 227ad6bc421SBarry Smith char type[256]; 228ad6bc421SBarry Smith 229ad6bc421SBarry Smith PetscFunctionBegin; 2304782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 231ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 232ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 233ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 234ad6bc421SBarry Smith 235060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2364782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2374782b174SLisandro Dalcin if (adapt->ops->load) { 2384782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 239ad6bc421SBarry Smith } 240ad6bc421SBarry Smith PetscFunctionReturn(0); 241ad6bc421SBarry Smith } 242ad6bc421SBarry Smith 24384df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 24484df9cb4SJed Brown { 24584df9cb4SJed Brown PetscErrorCode ierr; 2461c167fc2SEmil Constantinescu PetscBool iascii,isbinary,isnone,isglee; 24784df9cb4SJed Brown 24884df9cb4SJed Brown PetscFunctionBegin; 2494782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2504782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2514782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2524782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 253251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 254ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 25584df9cb4SJed Brown if (iascii) { 256dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 2571917a363SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr); 2581c167fc2SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTGLEE,&isglee);CHKERRQ(ierr); 2591917a363SLisandro Dalcin if (!isnone) { 260bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 2611917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 2621917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 2631917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 2641917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 265bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 266bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 2671c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," maximum solution absolute value to be ignored %g\n",(double)adapt->ignore_max);CHKERRQ(ierr); 2681c167fc2SEmil Constantinescu } 2691c167fc2SEmil Constantinescu if (isglee) { 2701c167fc2SEmil Constantinescu if (adapt->glee_use_local) { 2711c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE uses local error control\n");CHKERRQ(ierr); 2721c167fc2SEmil Constantinescu } else { 2731c167fc2SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE uses global error control\n");CHKERRQ(ierr); 2741c167fc2SEmil Constantinescu } 2751917a363SLisandro Dalcin } 27684df9cb4SJed Brown if (adapt->ops->view) { 27784df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 27884df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 27984df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 28084df9cb4SJed Brown } 281ad6bc421SBarry Smith } else if (isbinary) { 282ad6bc421SBarry Smith char type[256]; 283ad6bc421SBarry Smith 284ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 285ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 286ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 287bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 288f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 289f2c2a1b9SBarry Smith } 29084df9cb4SJed Brown PetscFunctionReturn(0); 29184df9cb4SJed Brown } 29284df9cb4SJed Brown 293099af0a3SLisandro Dalcin /*@ 294099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 295099af0a3SLisandro Dalcin 296099af0a3SLisandro Dalcin Collective on TS 297099af0a3SLisandro Dalcin 298099af0a3SLisandro Dalcin Input Parameter: 299099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 300099af0a3SLisandro Dalcin 301099af0a3SLisandro Dalcin Level: developer 302099af0a3SLisandro Dalcin 303099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 304099af0a3SLisandro Dalcin @*/ 305099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 306099af0a3SLisandro Dalcin { 307099af0a3SLisandro Dalcin PetscErrorCode ierr; 308099af0a3SLisandro Dalcin 309099af0a3SLisandro Dalcin PetscFunctionBegin; 310099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 311099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 312099af0a3SLisandro Dalcin PetscFunctionReturn(0); 313099af0a3SLisandro Dalcin } 314099af0a3SLisandro Dalcin 31584df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 31684df9cb4SJed Brown { 31784df9cb4SJed Brown PetscErrorCode ierr; 31884df9cb4SJed Brown 31984df9cb4SJed Brown PetscFunctionBegin; 32084df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 32184df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3224782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 323099af0a3SLisandro Dalcin 324099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 325099af0a3SLisandro Dalcin 32684df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3271c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 32884df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 32984df9cb4SJed Brown PetscFunctionReturn(0); 33084df9cb4SJed Brown } 33184df9cb4SJed Brown 3321c3436cfSJed Brown /*@ 3331c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3341c3436cfSJed Brown 3351c3436cfSJed Brown Collective on TSAdapt 3361c3436cfSJed Brown 3371c3436cfSJed Brown Input Arguments: 3381c3436cfSJed Brown + adapt - adaptive controller context 3391c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3401c3436cfSJed Brown 341bf997491SLisandro Dalcin Options Database Keys: 342ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 343bf997491SLisandro Dalcin 3441c3436cfSJed Brown Level: intermediate 3451c3436cfSJed Brown 3461c3436cfSJed Brown .seealso: TSAdaptChoose() 3471c3436cfSJed Brown @*/ 3481c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3491c3436cfSJed Brown { 3501c3436cfSJed Brown PetscErrorCode ierr; 3511c3436cfSJed Brown 3521c3436cfSJed Brown PetscFunctionBegin; 3534782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3544782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3551c3436cfSJed Brown if (flg) { 356ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3571c3436cfSJed Brown } else { 3581c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3591c3436cfSJed Brown } 3601c3436cfSJed Brown PetscFunctionReturn(0); 3611c3436cfSJed Brown } 3621c3436cfSJed Brown 3630873808bSJed Brown /*@C 364bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3650873808bSJed Brown 3660873808bSJed Brown Logically collective on TSAdapt 3670873808bSJed Brown 3680873808bSJed Brown Input Arguments: 3690873808bSJed Brown + adapt - adaptive controller context 3700873808bSJed Brown - func - stage check function 3710873808bSJed Brown 3720873808bSJed Brown Arguments of func: 3730873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3740873808bSJed Brown 3750873808bSJed Brown + adapt - adaptive controller context 3760873808bSJed Brown . ts - time stepping context 3770873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3780873808bSJed Brown 3790873808bSJed Brown Level: advanced 3800873808bSJed Brown 3810873808bSJed Brown .seealso: TSAdaptChoose() 3820873808bSJed Brown @*/ 383b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3840873808bSJed Brown { 3850873808bSJed Brown 3860873808bSJed Brown PetscFunctionBegin; 3870873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 38868ae941aSLisandro Dalcin adapt->checkstage = func; 3890873808bSJed Brown PetscFunctionReturn(0); 3900873808bSJed Brown } 3910873808bSJed Brown 3921c3436cfSJed Brown /*@ 393bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 394bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 395bf997491SLisandro Dalcin 3961917a363SLisandro Dalcin Logically collective on TSAdapt 397bf997491SLisandro Dalcin 398bf997491SLisandro Dalcin Input Arguments: 399bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 400bf997491SLisandro Dalcin - flag - whether to always accept steps 401bf997491SLisandro Dalcin 402bf997491SLisandro Dalcin Options Database Keys: 403ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 404bf997491SLisandro Dalcin 405bf997491SLisandro Dalcin Level: intermediate 406bf997491SLisandro Dalcin 407bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 408bf997491SLisandro Dalcin @*/ 409bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 410bf997491SLisandro Dalcin { 411bf997491SLisandro Dalcin PetscFunctionBegin; 412bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 413bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 414bf997491SLisandro Dalcin adapt->always_accept = flag; 415bf997491SLisandro Dalcin PetscFunctionReturn(0); 416bf997491SLisandro Dalcin } 417bf997491SLisandro Dalcin 418bf997491SLisandro Dalcin /*@ 4191917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4201c3436cfSJed Brown 4211917a363SLisandro Dalcin Logically collective on TSAdapt 4221917a363SLisandro Dalcin 4231917a363SLisandro Dalcin Input Arguments: 4241917a363SLisandro Dalcin + adapt - adaptive controller context 4251917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 426ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 4271917a363SLisandro Dalcin 4281917a363SLisandro Dalcin Options Database Keys: 429ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 430ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4311917a363SLisandro Dalcin 4321917a363SLisandro Dalcin Level: intermediate 4331917a363SLisandro Dalcin 4341917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4351917a363SLisandro Dalcin @*/ 4361917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4371917a363SLisandro Dalcin { 4381917a363SLisandro Dalcin PetscFunctionBegin; 4391917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4401917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4411917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4421917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4431917a363SLisandro 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); 4441917a363SLisandro 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); 4451917a363SLisandro 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); 4461917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4471917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4481917a363SLisandro Dalcin PetscFunctionReturn(0); 4491917a363SLisandro Dalcin } 4501917a363SLisandro Dalcin 4511917a363SLisandro Dalcin /*@ 4521917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4531917a363SLisandro Dalcin 4541917a363SLisandro Dalcin Not Collective 4551917a363SLisandro Dalcin 4561917a363SLisandro Dalcin Input Arguments: 4571917a363SLisandro Dalcin . adapt - adaptive controller context 4581917a363SLisandro Dalcin 4591917a363SLisandro Dalcin Ouput Arguments: 4601917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4611917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4621917a363SLisandro Dalcin 4631917a363SLisandro Dalcin Level: intermediate 4641917a363SLisandro Dalcin 4651917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4661917a363SLisandro Dalcin @*/ 4671917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4681917a363SLisandro Dalcin { 4691917a363SLisandro Dalcin PetscFunctionBegin; 4701917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4711917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4721917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4731917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4741917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4751917a363SLisandro Dalcin PetscFunctionReturn(0); 4761917a363SLisandro Dalcin } 4771917a363SLisandro Dalcin 4781917a363SLisandro Dalcin /*@ 479*76cddca1SEmil 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. 480*76cddca1SEmil Constantinescu 481*76cddca1SEmil Constantinescu Logically collective on TSAdapt 482*76cddca1SEmil Constantinescu 483*76cddca1SEmil Constantinescu Input Arguments: 484*76cddca1SEmil Constantinescu + adapt - adaptive controller context 485*76cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation 486*76cddca1SEmil Constantinescu 487*76cddca1SEmil Constantinescu Options Database Keys: 488*76cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold 489*76cddca1SEmil Constantinescu 490*76cddca1SEmil Constantinescu Level: intermediate 491*76cddca1SEmil Constantinescu 492*76cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptGetMaxIgnore(), TSAdaptChoose() 493*76cddca1SEmil Constantinescu @*/ 494*76cddca1SEmil Constantinescu PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore) 495*76cddca1SEmil Constantinescu { 496*76cddca1SEmil Constantinescu PetscFunctionBegin; 497*76cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 498*76cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt,max_ignore,2); 499*76cddca1SEmil Constantinescu adapt->ignore_max = max_ignore; 500*76cddca1SEmil Constantinescu PetscFunctionReturn(0); 501*76cddca1SEmil Constantinescu } 502*76cddca1SEmil Constantinescu 503*76cddca1SEmil Constantinescu /*@ 504*76cddca1SEmil 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). 505*76cddca1SEmil Constantinescu 506*76cddca1SEmil Constantinescu Not Collective 507*76cddca1SEmil Constantinescu 508*76cddca1SEmil Constantinescu Input Arguments: 509*76cddca1SEmil Constantinescu . adapt - adaptive controller context 510*76cddca1SEmil Constantinescu 511*76cddca1SEmil Constantinescu Ouput Arguments: 512*76cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation 513*76cddca1SEmil Constantinescu 514*76cddca1SEmil Constantinescu Level: intermediate 515*76cddca1SEmil Constantinescu 516*76cddca1SEmil Constantinescu .seealso: TSAdapt, TSAdaptSetMaxIgnore(), TSAdaptChoose() 517*76cddca1SEmil Constantinescu @*/ 518*76cddca1SEmil Constantinescu PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore) 519*76cddca1SEmil Constantinescu { 520*76cddca1SEmil Constantinescu PetscFunctionBegin; 521*76cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 522*76cddca1SEmil Constantinescu PetscValidRealPointer(max_ignore,2); 523*76cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max; 524*76cddca1SEmil Constantinescu PetscFunctionReturn(0); 525*76cddca1SEmil Constantinescu } 526*76cddca1SEmil Constantinescu 527*76cddca1SEmil Constantinescu 528*76cddca1SEmil Constantinescu /*@ 5291917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 5301917a363SLisandro Dalcin 5311917a363SLisandro Dalcin Logically collective on TSAdapt 5321917a363SLisandro Dalcin 5331917a363SLisandro Dalcin Input Arguments: 5341917a363SLisandro Dalcin + adapt - adaptive controller context 5351917a363SLisandro Dalcin . low - admissible decrease factor 536ec18b7bcSLisandro Dalcin - high - admissible increase factor 5371917a363SLisandro Dalcin 5381917a363SLisandro Dalcin Options Database Keys: 539ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 5401917a363SLisandro Dalcin 5411917a363SLisandro Dalcin Level: intermediate 5421917a363SLisandro Dalcin 5431917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptGetClip() 5441917a363SLisandro Dalcin @*/ 5451917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 5461917a363SLisandro Dalcin { 5471917a363SLisandro Dalcin PetscFunctionBegin; 5481917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5491917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 5501917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 5511917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 5521917a363SLisandro 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); 5531917a363SLisandro Dalcin if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high); 5541917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 5551917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 5561917a363SLisandro Dalcin PetscFunctionReturn(0); 5571917a363SLisandro Dalcin } 5581917a363SLisandro Dalcin 5591917a363SLisandro Dalcin /*@ 5601917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5611917a363SLisandro Dalcin 5621917a363SLisandro Dalcin Not Collective 5631917a363SLisandro Dalcin 5641917a363SLisandro Dalcin Input Arguments: 5651917a363SLisandro Dalcin . adapt - adaptive controller context 5661917a363SLisandro Dalcin 5671917a363SLisandro Dalcin Ouput Arguments: 5681917a363SLisandro Dalcin + low - optional, admissible decrease factor 5691917a363SLisandro Dalcin - high - optional, admissible increase factor 5701917a363SLisandro Dalcin 5711917a363SLisandro Dalcin Level: intermediate 5721917a363SLisandro Dalcin 5731917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptSetClip() 5741917a363SLisandro Dalcin @*/ 5751917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5761917a363SLisandro Dalcin { 5771917a363SLisandro Dalcin PetscFunctionBegin; 5781917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5791917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5801917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5811917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5821917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5831917a363SLisandro Dalcin PetscFunctionReturn(0); 5841917a363SLisandro Dalcin } 5851917a363SLisandro Dalcin 5861917a363SLisandro Dalcin /*@ 5871917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 5881917a363SLisandro Dalcin 5891917a363SLisandro Dalcin Logically collective on TSAdapt 5901c3436cfSJed Brown 5911c3436cfSJed Brown Input Arguments: 592552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5931c3436cfSJed Brown . hmin - minimum time step 5941c3436cfSJed Brown - hmax - maximum time step 5951c3436cfSJed Brown 5961c3436cfSJed Brown Options Database Keys: 597ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 598ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 5991c3436cfSJed Brown 6001c3436cfSJed Brown Level: intermediate 6011c3436cfSJed Brown 6021917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 6031c3436cfSJed Brown @*/ 6041c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 6051c3436cfSJed Brown { 6061c3436cfSJed Brown 6071c3436cfSJed Brown PetscFunctionBegin; 6084782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 609b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 610b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 611b1f5bccaSLisandro 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); 612b1f5bccaSLisandro 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); 613b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 614b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 615b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 616b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 617b1f5bccaSLisandro Dalcin if (hmax <= hmin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Maximum time step %g must geather than minimum time step %g",(double)hmax,(double)hmin); 6181917a363SLisandro Dalcin PetscFunctionReturn(0); 6191917a363SLisandro Dalcin } 6201917a363SLisandro Dalcin 6211917a363SLisandro Dalcin /*@ 6221917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 6231917a363SLisandro Dalcin 6241917a363SLisandro Dalcin Not Collective 6251917a363SLisandro Dalcin 6261917a363SLisandro Dalcin Input Arguments: 6271917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 6281917a363SLisandro Dalcin 6291917a363SLisandro Dalcin Output Arguments: 6301917a363SLisandro Dalcin + hmin - minimum time step 6311917a363SLisandro Dalcin - hmax - maximum time step 6321917a363SLisandro Dalcin 6331917a363SLisandro Dalcin Level: intermediate 6341917a363SLisandro Dalcin 6351917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 6361917a363SLisandro Dalcin @*/ 6371917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 6381917a363SLisandro Dalcin { 6391917a363SLisandro Dalcin 6401917a363SLisandro Dalcin PetscFunctionBegin; 6411917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6421917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 6431917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 6441917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 6451917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 6461c3436cfSJed Brown PetscFunctionReturn(0); 6471c3436cfSJed Brown } 6481c3436cfSJed Brown 649e55864a3SBarry Smith /* 65084df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 65184df9cb4SJed Brown 65284df9cb4SJed Brown Collective on TSAdapt 65384df9cb4SJed Brown 65484df9cb4SJed Brown Input Parameter: 65584df9cb4SJed Brown . adapt - the TSAdapt context 65684df9cb4SJed Brown 65784df9cb4SJed Brown Options Database Keys: 6581917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 659bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6601917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6611917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6621917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 66323de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 66423de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 66523de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 666de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 667de50f1caSBarry 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 66884df9cb4SJed Brown 66984df9cb4SJed Brown Level: advanced 67084df9cb4SJed Brown 67184df9cb4SJed Brown Notes: 67284df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 67384df9cb4SJed Brown 67423de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 67584df9cb4SJed Brown 6761917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 6771917a363SLisandro Dalcin TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 678e55864a3SBarry Smith */ 6794416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 68084df9cb4SJed Brown { 68184df9cb4SJed Brown PetscErrorCode ierr; 68284df9cb4SJed Brown char type[256] = TSADAPTBASIC; 6831917a363SLisandro Dalcin PetscReal safety,reject_safety,clip[2],hmin,hmax; 6841c3436cfSJed Brown PetscBool set,flg; 6851917a363SLisandro Dalcin PetscInt two; 68684df9cb4SJed Brown 68784df9cb4SJed Brown PetscFunctionBegin; 6884782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 68984df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6901566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 691e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 69223de1d84SBarry 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); 69384df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 69484df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 69584df9cb4SJed Brown } 6961917a363SLisandro Dalcin 697bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 698bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 6991917a363SLisandro Dalcin 7001917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 7011917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 7021917a363SLisandro 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); 7031917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 7041917a363SLisandro Dalcin 7051917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 7061917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 7071917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 7081917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 7091917a363SLisandro Dalcin 7101917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 7111917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 7121917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 713bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 7141917a363SLisandro Dalcin 715d580f011SEmil 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); 716d580f011SEmil 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); 717d580f011SEmil Constantinescu 7180298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);CHKERRQ(ierr); 7191917a363SLisandro Dalcin 7207619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 7217619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 7221917a363SLisandro Dalcin 723de50f1caSBarry 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); 724de50f1caSBarry Smith 7251917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 7261917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 7271917a363SLisandro Dalcin 728e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 72984df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 73084df9cb4SJed Brown PetscFunctionReturn(0); 73184df9cb4SJed Brown } 73284df9cb4SJed Brown 73384df9cb4SJed Brown /*@ 73484df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 73584df9cb4SJed Brown 7361917a363SLisandro Dalcin Logically collective on TSAdapt 73784df9cb4SJed Brown 73884df9cb4SJed Brown Input Argument: 73984df9cb4SJed Brown . adapt - adaptive controller 74084df9cb4SJed Brown 74184df9cb4SJed Brown Level: developer 74284df9cb4SJed Brown 74384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 74484df9cb4SJed Brown @*/ 74584df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 74684df9cb4SJed Brown { 74784df9cb4SJed Brown PetscErrorCode ierr; 74884df9cb4SJed Brown 74984df9cb4SJed Brown PetscFunctionBegin; 7504782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 75184df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 75284df9cb4SJed Brown PetscFunctionReturn(0); 75384df9cb4SJed Brown } 75484df9cb4SJed Brown 75584df9cb4SJed Brown /*@C 75684df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 75784df9cb4SJed Brown 7581917a363SLisandro Dalcin Logically collective on TSAdapt 75984df9cb4SJed Brown 76084df9cb4SJed Brown Input Arguments: 761552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 76284df9cb4SJed Brown . name - name of the candidate scheme to add 76384df9cb4SJed Brown . order - order of the candidate scheme 76484df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7658d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 76684df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 76784df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 76884df9cb4SJed Brown 76984df9cb4SJed Brown Note: 77084df9cb4SJed Brown This routine is not available in Fortran. 77184df9cb4SJed Brown 77284df9cb4SJed Brown Level: developer 77384df9cb4SJed Brown 77484df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 77584df9cb4SJed Brown @*/ 7768d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 77784df9cb4SJed Brown { 77884df9cb4SJed Brown PetscInt c; 77984df9cb4SJed Brown 78084df9cb4SJed Brown PetscFunctionBegin; 78184df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 782ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 78384df9cb4SJed Brown if (inuse) { 784ce94432eSBarry 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()"); 78584df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 78684df9cb4SJed Brown } 7871c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7881c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 789bbd56ea5SKarl Rupp 79084df9cb4SJed Brown adapt->candidates.name[c] = name; 79184df9cb4SJed Brown adapt->candidates.order[c] = order; 79284df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7938d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 79484df9cb4SJed Brown adapt->candidates.cost[c] = cost; 79584df9cb4SJed Brown adapt->candidates.n++; 79684df9cb4SJed Brown PetscFunctionReturn(0); 79784df9cb4SJed Brown } 79884df9cb4SJed Brown 7998d59e960SJed Brown /*@C 8008d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 8018d59e960SJed Brown 8028d59e960SJed Brown Not Collective 8038d59e960SJed Brown 8048d59e960SJed Brown Input Arguments: 8058d59e960SJed Brown . adapt - time step adaptivity context 8068d59e960SJed Brown 8078d59e960SJed Brown Output Arguments: 8088d59e960SJed Brown + n - number of candidate schemes, always at least 1 8098d59e960SJed Brown . order - the order of each candidate scheme 8108d59e960SJed Brown . stageorder - the stage order of each candidate scheme 8118d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 8128d59e960SJed Brown - cost - the relative cost of each scheme 8138d59e960SJed Brown 8148d59e960SJed Brown Level: developer 8158d59e960SJed Brown 8168d59e960SJed Brown Note: 8178d59e960SJed Brown The current scheme is always returned in the first slot 8188d59e960SJed Brown 8198d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 8208d59e960SJed Brown @*/ 8218d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 8228d59e960SJed Brown { 8238d59e960SJed Brown PetscFunctionBegin; 8248d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8258d59e960SJed Brown if (n) *n = adapt->candidates.n; 8268d59e960SJed Brown if (order) *order = adapt->candidates.order; 8278d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 8288d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 8298d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 8308d59e960SJed Brown PetscFunctionReturn(0); 8318d59e960SJed Brown } 8328d59e960SJed Brown 83384df9cb4SJed Brown /*@C 83484df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 83584df9cb4SJed Brown 8361917a363SLisandro Dalcin Collective on TSAdapt 83784df9cb4SJed Brown 83884df9cb4SJed Brown Input Arguments: 83984df9cb4SJed Brown + adapt - adaptive contoller 84084df9cb4SJed Brown - h - current step size 84184df9cb4SJed Brown 84284df9cb4SJed Brown Output Arguments: 8431566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 84484df9cb4SJed Brown . next_h - step size to use for the next step 84584df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 84684df9cb4SJed Brown 8471c3436cfSJed Brown Note: 8481c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 8491c3436cfSJed Brown being retried after an initial rejection. 8501c3436cfSJed Brown 85184df9cb4SJed Brown Level: developer 85284df9cb4SJed Brown 85384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 85484df9cb4SJed Brown @*/ 85584df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 85684df9cb4SJed Brown { 85784df9cb4SJed Brown PetscErrorCode ierr; 8581566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8591566a47fSLisandro Dalcin PetscInt scheme = 0; 8600b99f514SJed Brown PetscReal wlte = -1.0; 8615e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8625e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 86384df9cb4SJed Brown 86484df9cb4SJed Brown PetscFunctionBegin; 86584df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 86684df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8671566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 86884df9cb4SJed Brown PetscValidPointer(next_h,5); 86984df9cb4SJed Brown PetscValidIntPointer(accept,6); 8701566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8711566a47fSLisandro Dalcin 8721566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8731566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8741566a47fSLisandro Dalcin *next_h = h; 8751566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8761566a47fSLisandro Dalcin PetscFunctionReturn(0); 8771566a47fSLisandro Dalcin } 8781566a47fSLisandro Dalcin 8795e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 8801566a47fSLisandro 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); 8811566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8821566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8831566a47fSLisandro Dalcin 88449354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 88536b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 88636b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 88736b54a69SLisandro Dalcin PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 888fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 889fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 89036b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 891fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 89236b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 89349354f04SShri Abhyankar } 8941c3436cfSJed Brown 8951c3436cfSJed Brown if (adapt->monitor) { 8961566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8971c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8980b99f514SJed Brown if (wlte < 0) { 8990dea241dSBarry 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); 9000b99f514SJed Brown } else { 9010dea241dSBarry 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); 9020b99f514SJed Brown } 9031c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9041c3436cfSJed Brown } 90584df9cb4SJed Brown PetscFunctionReturn(0); 90684df9cb4SJed Brown } 90784df9cb4SJed Brown 90897335746SJed Brown /*@ 909de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 910de50f1caSBarry Smith before increasing the time step. 911de50f1caSBarry Smith 912de50f1caSBarry Smith Logicially Collective on TSAdapt 913de50f1caSBarry Smith 914de50f1caSBarry Smith Input Arguments: 915de50f1caSBarry Smith + adapt - adaptive controller context 916de50f1caSBarry Smith - cnt - the number of timesteps 917de50f1caSBarry Smith 918de50f1caSBarry Smith Options Database Key: 919de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 920de50f1caSBarry Smith 921de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 922de50f1caSBarry Smith The successful use of this option is problem dependent 923de50f1caSBarry Smith 924de50f1caSBarry Smith Developer Note: there is no theory to support this option 925de50f1caSBarry Smith 926de50f1caSBarry Smith Level: advanced 927de50f1caSBarry Smith 928de50f1caSBarry Smith .seealso: 929de50f1caSBarry Smith @*/ 930de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 931de50f1caSBarry Smith { 932de50f1caSBarry Smith PetscFunctionBegin; 933de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 934de50f1caSBarry Smith PetscFunctionReturn(0); 935de50f1caSBarry Smith } 936de50f1caSBarry Smith 937de50f1caSBarry Smith 938de50f1caSBarry Smith /*@ 93997335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 94097335746SJed Brown 9411917a363SLisandro Dalcin Collective on TSAdapt 94297335746SJed Brown 94397335746SJed Brown Input Arguments: 94497335746SJed Brown + adapt - adaptive controller context 945b295832fSPierre Barbier de Reuille . ts - time stepper 946b295832fSPierre Barbier de Reuille . t - Current simulation time 947b295832fSPierre Barbier de Reuille - Y - Current solution vector 94897335746SJed Brown 94997335746SJed Brown Output Arguments: 95097335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 95197335746SJed Brown 95297335746SJed Brown Level: developer 95397335746SJed Brown 95497335746SJed Brown .seealso: 95597335746SJed Brown @*/ 956b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 95797335746SJed Brown { 95897335746SJed Brown PetscErrorCode ierr; 9591566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 96097335746SJed Brown 96197335746SJed Brown PetscFunctionBegin; 9624782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 9634782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 9644782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 9651566a47fSLisandro Dalcin 9661566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 96797335746SJed Brown if (snesreason < 0) { 96897335746SJed Brown *accept = PETSC_FALSE; 9696de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 97097335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 97148c19aefSLisandro 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); 97297335746SJed Brown if (adapt->monitor) { 97397335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9740dea241dSBarry 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); 97597335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 97697335746SJed Brown } 977cb9d8021SPierre Barbier de Reuille } 978cb9d8021SPierre Barbier de Reuille } else { 9791566a47fSLisandro Dalcin *accept = PETSC_TRUE; 980b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 981cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 982b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 983cb9d8021SPierre Barbier de Reuille } 984cb9d8021SPierre Barbier de Reuille } 985cb9d8021SPierre Barbier de Reuille 9861566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 9871566a47fSLisandro Dalcin PetscReal dt,new_dt; 9881566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 989cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 99097335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 991de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 99297335746SJed Brown if (adapt->monitor) { 99397335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 994e6d0a238SBarry 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); 99597335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 99697335746SJed Brown } 99797335746SJed Brown } 99897335746SJed Brown PetscFunctionReturn(0); 99997335746SJed Brown } 100097335746SJed Brown 100184df9cb4SJed Brown /*@ 100284df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 100384df9cb4SJed Brown 100484df9cb4SJed Brown Collective on MPI_Comm 100584df9cb4SJed Brown 100684df9cb4SJed Brown Input Parameter: 100784df9cb4SJed Brown . comm - The communicator 100884df9cb4SJed Brown 100984df9cb4SJed Brown Output Parameter: 101084df9cb4SJed Brown . adapt - new TSAdapt object 101184df9cb4SJed Brown 101284df9cb4SJed Brown Level: developer 101384df9cb4SJed Brown 101484df9cb4SJed Brown Notes: 101584df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 101684df9cb4SJed Brown 101784df9cb4SJed Brown .keywords: TSAdapt, create 1018552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 101984df9cb4SJed Brown @*/ 102084df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 102184df9cb4SJed Brown { 102284df9cb4SJed Brown PetscErrorCode ierr; 102384df9cb4SJed Brown TSAdapt adapt; 102484df9cb4SJed Brown 102584df9cb4SJed Brown PetscFunctionBegin; 10263b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 10273b3bcf4cSLisandro Dalcin *inadapt = NULL; 10283b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 10293b3bcf4cSLisandro Dalcin 103073107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 10311c3436cfSJed Brown 1032bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 10331917a363SLisandro Dalcin adapt->safety = 0.9; 10341917a363SLisandro Dalcin adapt->reject_safety = 0.5; 10351917a363SLisandro Dalcin adapt->clip[0] = 0.1; 10361917a363SLisandro Dalcin adapt->clip[1] = 10.; 10371c3436cfSJed Brown adapt->dt_min = 1e-20; 10381917a363SLisandro Dalcin adapt->dt_max = 1e+20; 10391c167fc2SEmil Constantinescu adapt->ignore_max = -1.0; 1040d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 104197335746SJed Brown adapt->scale_solve_failed = 0.25; 1042fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1043fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1044fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1045fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 10467619abb3SShri adapt->wnormtype = NORM_2; 1047de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 10481917a363SLisandro Dalcin 104984df9cb4SJed Brown *inadapt = adapt; 105084df9cb4SJed Brown PetscFunctionReturn(0); 105184df9cb4SJed Brown } 1052