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 10484df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 1051917a363SLisandro Dalcin TSAdaptCreate() when using static libraries. 10684df9cb4SJed Brown 10784df9cb4SJed Brown Level: developer 10884df9cb4SJed Brown 10984df9cb4SJed Brown .keywords: TSAdapt, initialize, package 11084df9cb4SJed Brown .seealso: PetscInitialize() 11184df9cb4SJed Brown @*/ 112607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11384df9cb4SJed Brown { 11484df9cb4SJed Brown PetscErrorCode ierr; 11584df9cb4SJed Brown 11684df9cb4SJed Brown PetscFunctionBegin; 11784df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11884df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11984df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 120607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 12184df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 12284df9cb4SJed Brown PetscFunctionReturn(0); 12384df9cb4SJed Brown } 12484df9cb4SJed Brown 1257eef6055SBarry Smith /*@C 1267eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1277eef6055SBarry Smith 1287eef6055SBarry Smith Logicially Collective on TSAdapt 1297eef6055SBarry Smith 1307eef6055SBarry Smith Input Parameter: 131d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1327eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1337eef6055SBarry Smith 1347eef6055SBarry Smith Options Database: 135ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1367eef6055SBarry Smith 1377eef6055SBarry Smith Level: intermediate 1387eef6055SBarry Smith 1397eef6055SBarry Smith .keywords: TSAdapt, create 1407eef6055SBarry Smith 1417eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1427eef6055SBarry Smith @*/ 14319fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 14484df9cb4SJed Brown { 145ef749922SLisandro Dalcin PetscBool match; 14684df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 14784df9cb4SJed Brown 14884df9cb4SJed Brown PetscFunctionBegin; 1494782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 150b92453a8SLisandro Dalcin PetscValidCharPointer(type,2); 151ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 152ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1531c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 15484df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 155ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1564cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 15784df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15868ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15984df9cb4SJed Brown PetscFunctionReturn(0); 16084df9cb4SJed Brown } 16184df9cb4SJed Brown 162d0288e4fSLisandro Dalcin /*@C 163d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 164d0288e4fSLisandro Dalcin 165d0288e4fSLisandro Dalcin Not Collective 166d0288e4fSLisandro Dalcin 167d0288e4fSLisandro Dalcin Input Parameter: 168d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 169d0288e4fSLisandro Dalcin 170d0288e4fSLisandro Dalcin Output Parameter: 171d0288e4fSLisandro Dalcin . type - The name of TS adapter method 172d0288e4fSLisandro Dalcin 173d0288e4fSLisandro Dalcin Level: intermediate 174d0288e4fSLisandro Dalcin 175d0288e4fSLisandro Dalcin .keywords: TSAdapt, get, type 176d0288e4fSLisandro Dalcin .seealso TSAdaptSetType() 177d0288e4fSLisandro Dalcin @*/ 178d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 179d0288e4fSLisandro Dalcin { 180d0288e4fSLisandro Dalcin PetscFunctionBegin; 181d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 182d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 183d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 184d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 185d0288e4fSLisandro Dalcin } 186d0288e4fSLisandro Dalcin 18784df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 18884df9cb4SJed Brown { 18984df9cb4SJed Brown PetscErrorCode ierr; 19084df9cb4SJed Brown 19184df9cb4SJed Brown PetscFunctionBegin; 1924782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 19384df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 19484df9cb4SJed Brown PetscFunctionReturn(0); 19584df9cb4SJed Brown } 19684df9cb4SJed Brown 197ad6bc421SBarry Smith /*@C 198ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 199ad6bc421SBarry Smith 200ad6bc421SBarry Smith Collective on PetscViewer 201ad6bc421SBarry Smith 202ad6bc421SBarry Smith Input Parameters: 203ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 204ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 205ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 206ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 207ad6bc421SBarry Smith 208ad6bc421SBarry Smith Level: intermediate 209ad6bc421SBarry Smith 210ad6bc421SBarry Smith Notes: 211ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 212ad6bc421SBarry Smith 213ad6bc421SBarry Smith Notes for advanced users: 214ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 215ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 216ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 217ad6bc421SBarry Smith format is 218ad6bc421SBarry Smith .vb 219ad6bc421SBarry Smith has not yet been determined 220ad6bc421SBarry Smith .ve 221ad6bc421SBarry Smith 222ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 223ad6bc421SBarry Smith @*/ 2244782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 225ad6bc421SBarry Smith { 226ad6bc421SBarry Smith PetscErrorCode ierr; 227ad6bc421SBarry Smith PetscBool isbinary; 228ad6bc421SBarry Smith char type[256]; 229ad6bc421SBarry Smith 230ad6bc421SBarry Smith PetscFunctionBegin; 2314782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 232ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 233ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 234ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 235ad6bc421SBarry Smith 236060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2374782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2384782b174SLisandro Dalcin if (adapt->ops->load) { 2394782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 240ad6bc421SBarry Smith } 241ad6bc421SBarry Smith PetscFunctionReturn(0); 242ad6bc421SBarry Smith } 243ad6bc421SBarry Smith 24484df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 24584df9cb4SJed Brown { 24684df9cb4SJed Brown PetscErrorCode ierr; 2471917a363SLisandro Dalcin PetscBool iascii,isbinary,isnone; 24884df9cb4SJed Brown 24984df9cb4SJed Brown PetscFunctionBegin; 2504782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2514782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2524782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2534782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 254251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 255ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 25684df9cb4SJed Brown if (iascii) { 257dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 2581917a363SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);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); 2671917a363SLisandro Dalcin } 26884df9cb4SJed Brown if (adapt->ops->view) { 26984df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 27084df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 27184df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27284df9cb4SJed Brown } 273ad6bc421SBarry Smith } else if (isbinary) { 274ad6bc421SBarry Smith char type[256]; 275ad6bc421SBarry Smith 276ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 277ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 278ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 279bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 280f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 281f2c2a1b9SBarry Smith } 28284df9cb4SJed Brown PetscFunctionReturn(0); 28384df9cb4SJed Brown } 28484df9cb4SJed Brown 285099af0a3SLisandro Dalcin /*@ 286099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 287099af0a3SLisandro Dalcin 288099af0a3SLisandro Dalcin Collective on TS 289099af0a3SLisandro Dalcin 290099af0a3SLisandro Dalcin Input Parameter: 291099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 292099af0a3SLisandro Dalcin 293099af0a3SLisandro Dalcin Level: developer 294099af0a3SLisandro Dalcin 295099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 296099af0a3SLisandro Dalcin @*/ 297099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 298099af0a3SLisandro Dalcin { 299099af0a3SLisandro Dalcin PetscErrorCode ierr; 300099af0a3SLisandro Dalcin 301099af0a3SLisandro Dalcin PetscFunctionBegin; 302099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 303099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 304099af0a3SLisandro Dalcin PetscFunctionReturn(0); 305099af0a3SLisandro Dalcin } 306099af0a3SLisandro Dalcin 30784df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 30884df9cb4SJed Brown { 30984df9cb4SJed Brown PetscErrorCode ierr; 31084df9cb4SJed Brown 31184df9cb4SJed Brown PetscFunctionBegin; 31284df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 31384df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3144782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 315099af0a3SLisandro Dalcin 316099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 317099af0a3SLisandro Dalcin 31884df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3191c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 32084df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 32184df9cb4SJed Brown PetscFunctionReturn(0); 32284df9cb4SJed Brown } 32384df9cb4SJed Brown 3241c3436cfSJed Brown /*@ 3251c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3261c3436cfSJed Brown 3271c3436cfSJed Brown Collective on TSAdapt 3281c3436cfSJed Brown 3291c3436cfSJed Brown Input Arguments: 3301c3436cfSJed Brown + adapt - adaptive controller context 3311c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3321c3436cfSJed Brown 333bf997491SLisandro Dalcin Options Database Keys: 334ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 335bf997491SLisandro Dalcin 3361c3436cfSJed Brown Level: intermediate 3371c3436cfSJed Brown 3381c3436cfSJed Brown .seealso: TSAdaptChoose() 3391c3436cfSJed Brown @*/ 3401c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3411c3436cfSJed Brown { 3421c3436cfSJed Brown PetscErrorCode ierr; 3431c3436cfSJed Brown 3441c3436cfSJed Brown PetscFunctionBegin; 3454782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3464782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3471c3436cfSJed Brown if (flg) { 348ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3491c3436cfSJed Brown } else { 3501c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3511c3436cfSJed Brown } 3521c3436cfSJed Brown PetscFunctionReturn(0); 3531c3436cfSJed Brown } 3541c3436cfSJed Brown 3550873808bSJed Brown /*@C 356bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3570873808bSJed Brown 3580873808bSJed Brown Logically collective on TSAdapt 3590873808bSJed Brown 3600873808bSJed Brown Input Arguments: 3610873808bSJed Brown + adapt - adaptive controller context 3620873808bSJed Brown - func - stage check function 3630873808bSJed Brown 3640873808bSJed Brown Arguments of func: 3650873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3660873808bSJed Brown 3670873808bSJed Brown + adapt - adaptive controller context 3680873808bSJed Brown . ts - time stepping context 3690873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3700873808bSJed Brown 3710873808bSJed Brown Level: advanced 3720873808bSJed Brown 3730873808bSJed Brown .seealso: TSAdaptChoose() 3740873808bSJed Brown @*/ 375b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3760873808bSJed Brown { 3770873808bSJed Brown 3780873808bSJed Brown PetscFunctionBegin; 3790873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 38068ae941aSLisandro Dalcin adapt->checkstage = func; 3810873808bSJed Brown PetscFunctionReturn(0); 3820873808bSJed Brown } 3830873808bSJed Brown 3841c3436cfSJed Brown /*@ 385bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 386bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 387bf997491SLisandro Dalcin 3881917a363SLisandro Dalcin Logically collective on TSAdapt 389bf997491SLisandro Dalcin 390bf997491SLisandro Dalcin Input Arguments: 391bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 392bf997491SLisandro Dalcin - flag - whether to always accept steps 393bf997491SLisandro Dalcin 394bf997491SLisandro Dalcin Options Database Keys: 395ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 396bf997491SLisandro Dalcin 397bf997491SLisandro Dalcin Level: intermediate 398bf997491SLisandro Dalcin 399bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 400bf997491SLisandro Dalcin @*/ 401bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 402bf997491SLisandro Dalcin { 403bf997491SLisandro Dalcin PetscFunctionBegin; 404bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 405bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 406bf997491SLisandro Dalcin adapt->always_accept = flag; 407bf997491SLisandro Dalcin PetscFunctionReturn(0); 408bf997491SLisandro Dalcin } 409bf997491SLisandro Dalcin 410bf997491SLisandro Dalcin /*@ 4111917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4121c3436cfSJed Brown 4131917a363SLisandro Dalcin Logically collective on TSAdapt 4141917a363SLisandro Dalcin 4151917a363SLisandro Dalcin Input Arguments: 4161917a363SLisandro Dalcin + adapt - adaptive controller context 4171917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 418ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 4191917a363SLisandro Dalcin 4201917a363SLisandro Dalcin Options Database Keys: 421ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 422ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4231917a363SLisandro Dalcin 4241917a363SLisandro Dalcin Level: intermediate 4251917a363SLisandro Dalcin 4261917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4271917a363SLisandro Dalcin @*/ 4281917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4291917a363SLisandro Dalcin { 4301917a363SLisandro Dalcin PetscFunctionBegin; 4311917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4321917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4331917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4341917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4351917a363SLisandro 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); 4361917a363SLisandro 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); 4371917a363SLisandro 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); 4381917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4391917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4401917a363SLisandro Dalcin PetscFunctionReturn(0); 4411917a363SLisandro Dalcin } 4421917a363SLisandro Dalcin 4431917a363SLisandro Dalcin /*@ 4441917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4451917a363SLisandro Dalcin 4461917a363SLisandro Dalcin Not Collective 4471917a363SLisandro Dalcin 4481917a363SLisandro Dalcin Input Arguments: 4491917a363SLisandro Dalcin . adapt - adaptive controller context 4501917a363SLisandro Dalcin 4511917a363SLisandro Dalcin Ouput Arguments: 4521917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4531917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4541917a363SLisandro Dalcin 4551917a363SLisandro Dalcin Level: intermediate 4561917a363SLisandro Dalcin 4571917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4581917a363SLisandro Dalcin @*/ 4591917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4601917a363SLisandro Dalcin { 4611917a363SLisandro Dalcin PetscFunctionBegin; 4621917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4631917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4641917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4651917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4661917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4671917a363SLisandro Dalcin PetscFunctionReturn(0); 4681917a363SLisandro Dalcin } 4691917a363SLisandro Dalcin 4701917a363SLisandro Dalcin /*@ 4711917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4721917a363SLisandro Dalcin 4731917a363SLisandro Dalcin Logically collective on TSAdapt 4741917a363SLisandro Dalcin 4751917a363SLisandro Dalcin Input Arguments: 4761917a363SLisandro Dalcin + adapt - adaptive controller context 4771917a363SLisandro Dalcin . low - admissible decrease factor 478ec18b7bcSLisandro Dalcin - high - admissible increase factor 4791917a363SLisandro Dalcin 4801917a363SLisandro Dalcin Options Database Keys: 481ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 4821917a363SLisandro Dalcin 4831917a363SLisandro Dalcin Level: intermediate 4841917a363SLisandro Dalcin 4851917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptGetClip() 4861917a363SLisandro Dalcin @*/ 4871917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 4881917a363SLisandro Dalcin { 4891917a363SLisandro Dalcin PetscFunctionBegin; 4901917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4911917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 4921917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 4931917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 4941917a363SLisandro 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); 4951917a363SLisandro 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); 4961917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 4971917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 4981917a363SLisandro Dalcin PetscFunctionReturn(0); 4991917a363SLisandro Dalcin } 5001917a363SLisandro Dalcin 5011917a363SLisandro Dalcin /*@ 5021917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5031917a363SLisandro Dalcin 5041917a363SLisandro Dalcin Not Collective 5051917a363SLisandro Dalcin 5061917a363SLisandro Dalcin Input Arguments: 5071917a363SLisandro Dalcin . adapt - adaptive controller context 5081917a363SLisandro Dalcin 5091917a363SLisandro Dalcin Ouput Arguments: 5101917a363SLisandro Dalcin + low - optional, admissible decrease factor 5111917a363SLisandro Dalcin - high - optional, admissible increase factor 5121917a363SLisandro Dalcin 5131917a363SLisandro Dalcin Level: intermediate 5141917a363SLisandro Dalcin 5151917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptSetClip() 5161917a363SLisandro Dalcin @*/ 5171917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5181917a363SLisandro Dalcin { 5191917a363SLisandro Dalcin PetscFunctionBegin; 5201917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5211917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5221917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5231917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5241917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5251917a363SLisandro Dalcin PetscFunctionReturn(0); 5261917a363SLisandro Dalcin } 5271917a363SLisandro Dalcin 5281917a363SLisandro Dalcin /*@ 5291917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 5301917a363SLisandro Dalcin 5311917a363SLisandro Dalcin Logically collective on TSAdapt 5321c3436cfSJed Brown 5331c3436cfSJed Brown Input Arguments: 534552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5351c3436cfSJed Brown . hmin - minimum time step 5361c3436cfSJed Brown - hmax - maximum time step 5371c3436cfSJed Brown 5381c3436cfSJed Brown Options Database Keys: 539ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 540ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 5411c3436cfSJed Brown 5421c3436cfSJed Brown Level: intermediate 5431c3436cfSJed Brown 5441917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 5451c3436cfSJed Brown @*/ 5461c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 5471c3436cfSJed Brown { 5481c3436cfSJed Brown 5491c3436cfSJed Brown PetscFunctionBegin; 5504782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 551b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 552b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 553b1f5bccaSLisandro 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); 554b1f5bccaSLisandro 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); 555b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 556b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 557b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 558b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 559b1f5bccaSLisandro 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); 5601917a363SLisandro Dalcin PetscFunctionReturn(0); 5611917a363SLisandro Dalcin } 5621917a363SLisandro Dalcin 5631917a363SLisandro Dalcin /*@ 5641917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 5651917a363SLisandro Dalcin 5661917a363SLisandro Dalcin Not Collective 5671917a363SLisandro Dalcin 5681917a363SLisandro Dalcin Input Arguments: 5691917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5701917a363SLisandro Dalcin 5711917a363SLisandro Dalcin Output Arguments: 5721917a363SLisandro Dalcin + hmin - minimum time step 5731917a363SLisandro Dalcin - hmax - maximum time step 5741917a363SLisandro Dalcin 5751917a363SLisandro Dalcin Level: intermediate 5761917a363SLisandro Dalcin 5771917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 5781917a363SLisandro Dalcin @*/ 5791917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 5801917a363SLisandro Dalcin { 5811917a363SLisandro Dalcin 5821917a363SLisandro Dalcin PetscFunctionBegin; 5831917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5841917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 5851917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 5861917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 5871917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 5881c3436cfSJed Brown PetscFunctionReturn(0); 5891c3436cfSJed Brown } 5901c3436cfSJed Brown 591e55864a3SBarry Smith /* 59284df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 59384df9cb4SJed Brown 59484df9cb4SJed Brown Collective on TSAdapt 59584df9cb4SJed Brown 59684df9cb4SJed Brown Input Parameter: 59784df9cb4SJed Brown . adapt - the TSAdapt context 59884df9cb4SJed Brown 59984df9cb4SJed Brown Options Database Keys: 6001917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 601bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6021917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6031917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6041917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 60523de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 60623de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 60723de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 608*de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 609*de50f1caSBarry 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 61084df9cb4SJed Brown 61184df9cb4SJed Brown Level: advanced 61284df9cb4SJed Brown 61384df9cb4SJed Brown Notes: 61484df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 61584df9cb4SJed Brown 61623de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 61784df9cb4SJed Brown 6181917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 6191917a363SLisandro Dalcin TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 620e55864a3SBarry Smith */ 6214416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 62284df9cb4SJed Brown { 62384df9cb4SJed Brown PetscErrorCode ierr; 62484df9cb4SJed Brown char type[256] = TSADAPTBASIC; 6251917a363SLisandro Dalcin PetscReal safety,reject_safety,clip[2],hmin,hmax; 6261c3436cfSJed Brown PetscBool set,flg; 6271917a363SLisandro Dalcin PetscInt two; 62884df9cb4SJed Brown 62984df9cb4SJed Brown PetscFunctionBegin; 6304782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 63184df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6321566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 633e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 63423de1d84SBarry 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); 63584df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 63684df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 63784df9cb4SJed Brown } 6381917a363SLisandro Dalcin 639bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 640bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 6411917a363SLisandro Dalcin 6421917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 6431917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 6441917a363SLisandro 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); 6451917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 6461917a363SLisandro Dalcin 6471917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 6481917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 6491917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 6501917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 6511917a363SLisandro Dalcin 6521917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 6531917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 6541917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 655bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 6561917a363SLisandro Dalcin 6570298fd71SBarry 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); 6581917a363SLisandro Dalcin 6597619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 6607619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 6611917a363SLisandro Dalcin 662*de50f1caSBarry 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); 663*de50f1caSBarry Smith 6641917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 6651917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 6661917a363SLisandro Dalcin 667e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 66884df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 66984df9cb4SJed Brown PetscFunctionReturn(0); 67084df9cb4SJed Brown } 67184df9cb4SJed Brown 67284df9cb4SJed Brown /*@ 67384df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 67484df9cb4SJed Brown 6751917a363SLisandro Dalcin Logically collective on TSAdapt 67684df9cb4SJed Brown 67784df9cb4SJed Brown Input Argument: 67884df9cb4SJed Brown . adapt - adaptive controller 67984df9cb4SJed Brown 68084df9cb4SJed Brown Level: developer 68184df9cb4SJed Brown 68284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 68384df9cb4SJed Brown @*/ 68484df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 68584df9cb4SJed Brown { 68684df9cb4SJed Brown PetscErrorCode ierr; 68784df9cb4SJed Brown 68884df9cb4SJed Brown PetscFunctionBegin; 6894782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 69084df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 69184df9cb4SJed Brown PetscFunctionReturn(0); 69284df9cb4SJed Brown } 69384df9cb4SJed Brown 69484df9cb4SJed Brown /*@C 69584df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 69684df9cb4SJed Brown 6971917a363SLisandro Dalcin Logically collective on TSAdapt 69884df9cb4SJed Brown 69984df9cb4SJed Brown Input Arguments: 700552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 70184df9cb4SJed Brown . name - name of the candidate scheme to add 70284df9cb4SJed Brown . order - order of the candidate scheme 70384df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7048d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 70584df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 70684df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 70784df9cb4SJed Brown 70884df9cb4SJed Brown Note: 70984df9cb4SJed Brown This routine is not available in Fortran. 71084df9cb4SJed Brown 71184df9cb4SJed Brown Level: developer 71284df9cb4SJed Brown 71384df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 71484df9cb4SJed Brown @*/ 7158d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 71684df9cb4SJed Brown { 71784df9cb4SJed Brown PetscInt c; 71884df9cb4SJed Brown 71984df9cb4SJed Brown PetscFunctionBegin; 72084df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 721ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 72284df9cb4SJed Brown if (inuse) { 723ce94432eSBarry 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()"); 72484df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 72584df9cb4SJed Brown } 7261c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7271c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 728bbd56ea5SKarl Rupp 72984df9cb4SJed Brown adapt->candidates.name[c] = name; 73084df9cb4SJed Brown adapt->candidates.order[c] = order; 73184df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7328d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 73384df9cb4SJed Brown adapt->candidates.cost[c] = cost; 73484df9cb4SJed Brown adapt->candidates.n++; 73584df9cb4SJed Brown PetscFunctionReturn(0); 73684df9cb4SJed Brown } 73784df9cb4SJed Brown 7388d59e960SJed Brown /*@C 7398d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 7408d59e960SJed Brown 7418d59e960SJed Brown Not Collective 7428d59e960SJed Brown 7438d59e960SJed Brown Input Arguments: 7448d59e960SJed Brown . adapt - time step adaptivity context 7458d59e960SJed Brown 7468d59e960SJed Brown Output Arguments: 7478d59e960SJed Brown + n - number of candidate schemes, always at least 1 7488d59e960SJed Brown . order - the order of each candidate scheme 7498d59e960SJed Brown . stageorder - the stage order of each candidate scheme 7508d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 7518d59e960SJed Brown - cost - the relative cost of each scheme 7528d59e960SJed Brown 7538d59e960SJed Brown Level: developer 7548d59e960SJed Brown 7558d59e960SJed Brown Note: 7568d59e960SJed Brown The current scheme is always returned in the first slot 7578d59e960SJed Brown 7588d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 7598d59e960SJed Brown @*/ 7608d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 7618d59e960SJed Brown { 7628d59e960SJed Brown PetscFunctionBegin; 7638d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7648d59e960SJed Brown if (n) *n = adapt->candidates.n; 7658d59e960SJed Brown if (order) *order = adapt->candidates.order; 7668d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 7678d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 7688d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 7698d59e960SJed Brown PetscFunctionReturn(0); 7708d59e960SJed Brown } 7718d59e960SJed Brown 77284df9cb4SJed Brown /*@C 77384df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 77484df9cb4SJed Brown 7751917a363SLisandro Dalcin Collective on TSAdapt 77684df9cb4SJed Brown 77784df9cb4SJed Brown Input Arguments: 77884df9cb4SJed Brown + adapt - adaptive contoller 77984df9cb4SJed Brown - h - current step size 78084df9cb4SJed Brown 78184df9cb4SJed Brown Output Arguments: 7821566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 78384df9cb4SJed Brown . next_h - step size to use for the next step 78484df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 78584df9cb4SJed Brown 7861c3436cfSJed Brown Note: 7871c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 7881c3436cfSJed Brown being retried after an initial rejection. 7891c3436cfSJed Brown 79084df9cb4SJed Brown Level: developer 79184df9cb4SJed Brown 79284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 79384df9cb4SJed Brown @*/ 79484df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 79584df9cb4SJed Brown { 79684df9cb4SJed Brown PetscErrorCode ierr; 7971566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 7981566a47fSLisandro Dalcin PetscInt scheme = 0; 7990b99f514SJed Brown PetscReal wlte = -1.0; 8005e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8015e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 80284df9cb4SJed Brown 80384df9cb4SJed Brown PetscFunctionBegin; 80484df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 80584df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8061566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 80784df9cb4SJed Brown PetscValidPointer(next_h,5); 80884df9cb4SJed Brown PetscValidIntPointer(accept,6); 8091566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8101566a47fSLisandro Dalcin 8111566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8121566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8131566a47fSLisandro Dalcin *next_h = h; 8141566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8151566a47fSLisandro Dalcin PetscFunctionReturn(0); 8161566a47fSLisandro Dalcin } 8171566a47fSLisandro Dalcin 8185e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 8191566a47fSLisandro 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); 8201566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8211566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8221566a47fSLisandro Dalcin 82349354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 82436b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 82536b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 82636b54a69SLisandro Dalcin PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 827fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 828fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 82936b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 830fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 83136b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 83249354f04SShri Abhyankar } 8331c3436cfSJed Brown 8341c3436cfSJed Brown if (adapt->monitor) { 8351566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8361c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8370b99f514SJed Brown if (wlte < 0) { 8380dea241dSBarry 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); 8390b99f514SJed Brown } else { 8400dea241dSBarry 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); 8410b99f514SJed Brown } 8421c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8431c3436cfSJed Brown } 84484df9cb4SJed Brown PetscFunctionReturn(0); 84584df9cb4SJed Brown } 84684df9cb4SJed Brown 84797335746SJed Brown /*@ 848*de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 849*de50f1caSBarry Smith before increasing the time step. 850*de50f1caSBarry Smith 851*de50f1caSBarry Smith Logicially Collective on TSAdapt 852*de50f1caSBarry Smith 853*de50f1caSBarry Smith Input Arguments: 854*de50f1caSBarry Smith + adapt - adaptive controller context 855*de50f1caSBarry Smith - cnt - the number of timesteps 856*de50f1caSBarry Smith 857*de50f1caSBarry Smith Options Database Key: 858*de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 859*de50f1caSBarry Smith 860*de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 861*de50f1caSBarry Smith The successful use of this option is problem dependent 862*de50f1caSBarry Smith 863*de50f1caSBarry Smith Developer Note: there is no theory to support this option 864*de50f1caSBarry Smith 865*de50f1caSBarry Smith Level: advanced 866*de50f1caSBarry Smith 867*de50f1caSBarry Smith .seealso: 868*de50f1caSBarry Smith @*/ 869*de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 870*de50f1caSBarry Smith { 871*de50f1caSBarry Smith PetscFunctionBegin; 872*de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 873*de50f1caSBarry Smith PetscFunctionReturn(0); 874*de50f1caSBarry Smith } 875*de50f1caSBarry Smith 876*de50f1caSBarry Smith 877*de50f1caSBarry Smith /*@ 87897335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 87997335746SJed Brown 8801917a363SLisandro Dalcin Collective on TSAdapt 88197335746SJed Brown 88297335746SJed Brown Input Arguments: 88397335746SJed Brown + adapt - adaptive controller context 884b295832fSPierre Barbier de Reuille . ts - time stepper 885b295832fSPierre Barbier de Reuille . t - Current simulation time 886b295832fSPierre Barbier de Reuille - Y - Current solution vector 88797335746SJed Brown 88897335746SJed Brown Output Arguments: 88997335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 89097335746SJed Brown 89197335746SJed Brown Level: developer 89297335746SJed Brown 89397335746SJed Brown .seealso: 89497335746SJed Brown @*/ 895b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 89697335746SJed Brown { 89797335746SJed Brown PetscErrorCode ierr; 8981566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 89997335746SJed Brown 90097335746SJed Brown PetscFunctionBegin; 9014782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 9024782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 9034782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 9041566a47fSLisandro Dalcin 9051566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 90697335746SJed Brown if (snesreason < 0) { 90797335746SJed Brown *accept = PETSC_FALSE; 9086de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 90997335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 91048c19aefSLisandro 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); 91197335746SJed Brown if (adapt->monitor) { 91297335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9130dea241dSBarry 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); 91497335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 91597335746SJed Brown } 916cb9d8021SPierre Barbier de Reuille } 917cb9d8021SPierre Barbier de Reuille } else { 9181566a47fSLisandro Dalcin *accept = PETSC_TRUE; 919b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 920cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 921b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 922cb9d8021SPierre Barbier de Reuille } 923cb9d8021SPierre Barbier de Reuille } 924cb9d8021SPierre Barbier de Reuille 9251566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 9261566a47fSLisandro Dalcin PetscReal dt,new_dt; 9271566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 928cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 92997335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 930*de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 93197335746SJed Brown if (adapt->monitor) { 93297335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 933e6d0a238SBarry 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); 93497335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 93597335746SJed Brown } 93697335746SJed Brown } 93797335746SJed Brown PetscFunctionReturn(0); 93897335746SJed Brown } 93997335746SJed Brown 94084df9cb4SJed Brown /*@ 94184df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 94284df9cb4SJed Brown 94384df9cb4SJed Brown Collective on MPI_Comm 94484df9cb4SJed Brown 94584df9cb4SJed Brown Input Parameter: 94684df9cb4SJed Brown . comm - The communicator 94784df9cb4SJed Brown 94884df9cb4SJed Brown Output Parameter: 94984df9cb4SJed Brown . adapt - new TSAdapt object 95084df9cb4SJed Brown 95184df9cb4SJed Brown Level: developer 95284df9cb4SJed Brown 95384df9cb4SJed Brown Notes: 95484df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 95584df9cb4SJed Brown 95684df9cb4SJed Brown .keywords: TSAdapt, create 957552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 95884df9cb4SJed Brown @*/ 95984df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 96084df9cb4SJed Brown { 96184df9cb4SJed Brown PetscErrorCode ierr; 96284df9cb4SJed Brown TSAdapt adapt; 96384df9cb4SJed Brown 96484df9cb4SJed Brown PetscFunctionBegin; 9653b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 9663b3bcf4cSLisandro Dalcin *inadapt = NULL; 9673b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 9683b3bcf4cSLisandro Dalcin 96973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 9701c3436cfSJed Brown 971bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 9721917a363SLisandro Dalcin adapt->safety = 0.9; 9731917a363SLisandro Dalcin adapt->reject_safety = 0.5; 9741917a363SLisandro Dalcin adapt->clip[0] = 0.1; 9751917a363SLisandro Dalcin adapt->clip[1] = 10.; 9761c3436cfSJed Brown adapt->dt_min = 1e-20; 9771917a363SLisandro Dalcin adapt->dt_max = 1e+20; 97897335746SJed Brown adapt->scale_solve_failed = 0.25; 979fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 980fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 981fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 982fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 9837619abb3SShri adapt->wnormtype = NORM_2; 984*de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 9851917a363SLisandro Dalcin 98684df9cb4SJed Brown *inadapt = adapt; 98784df9cb4SJed Brown PetscFunctionReturn(0); 98884df9cb4SJed Brown } 989