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; 2461917a363SLisandro Dalcin PetscBool iascii,isbinary,isnone; 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); 2581917a363SLisandro Dalcin if (!isnone) { 259bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 2601917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 2611917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 2621917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 2631917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 264bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 265bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 2661917a363SLisandro Dalcin } 26784df9cb4SJed Brown if (adapt->ops->view) { 26884df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 26984df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 27084df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27184df9cb4SJed Brown } 272ad6bc421SBarry Smith } else if (isbinary) { 273ad6bc421SBarry Smith char type[256]; 274ad6bc421SBarry Smith 275ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 276ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 277ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 278bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 279f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 280f2c2a1b9SBarry Smith } 28184df9cb4SJed Brown PetscFunctionReturn(0); 28284df9cb4SJed Brown } 28384df9cb4SJed Brown 284099af0a3SLisandro Dalcin /*@ 285099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 286099af0a3SLisandro Dalcin 287099af0a3SLisandro Dalcin Collective on TS 288099af0a3SLisandro Dalcin 289099af0a3SLisandro Dalcin Input Parameter: 290099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 291099af0a3SLisandro Dalcin 292099af0a3SLisandro Dalcin Level: developer 293099af0a3SLisandro Dalcin 294099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 295099af0a3SLisandro Dalcin @*/ 296099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 297099af0a3SLisandro Dalcin { 298099af0a3SLisandro Dalcin PetscErrorCode ierr; 299099af0a3SLisandro Dalcin 300099af0a3SLisandro Dalcin PetscFunctionBegin; 301099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 302099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 303099af0a3SLisandro Dalcin PetscFunctionReturn(0); 304099af0a3SLisandro Dalcin } 305099af0a3SLisandro Dalcin 30684df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 30784df9cb4SJed Brown { 30884df9cb4SJed Brown PetscErrorCode ierr; 30984df9cb4SJed Brown 31084df9cb4SJed Brown PetscFunctionBegin; 31184df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 31284df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3134782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 314099af0a3SLisandro Dalcin 315099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 316099af0a3SLisandro Dalcin 31784df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3181c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 31984df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 32084df9cb4SJed Brown PetscFunctionReturn(0); 32184df9cb4SJed Brown } 32284df9cb4SJed Brown 3231c3436cfSJed Brown /*@ 3241c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3251c3436cfSJed Brown 3261c3436cfSJed Brown Collective on TSAdapt 3271c3436cfSJed Brown 3281c3436cfSJed Brown Input Arguments: 3291c3436cfSJed Brown + adapt - adaptive controller context 3301c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3311c3436cfSJed Brown 332bf997491SLisandro Dalcin Options Database Keys: 333ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 334bf997491SLisandro Dalcin 3351c3436cfSJed Brown Level: intermediate 3361c3436cfSJed Brown 3371c3436cfSJed Brown .seealso: TSAdaptChoose() 3381c3436cfSJed Brown @*/ 3391c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3401c3436cfSJed Brown { 3411c3436cfSJed Brown PetscErrorCode ierr; 3421c3436cfSJed Brown 3431c3436cfSJed Brown PetscFunctionBegin; 3444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3454782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3461c3436cfSJed Brown if (flg) { 347ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3481c3436cfSJed Brown } else { 3491c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3501c3436cfSJed Brown } 3511c3436cfSJed Brown PetscFunctionReturn(0); 3521c3436cfSJed Brown } 3531c3436cfSJed Brown 3540873808bSJed Brown /*@C 355bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3560873808bSJed Brown 3570873808bSJed Brown Logically collective on TSAdapt 3580873808bSJed Brown 3590873808bSJed Brown Input Arguments: 3600873808bSJed Brown + adapt - adaptive controller context 3610873808bSJed Brown - func - stage check function 3620873808bSJed Brown 3630873808bSJed Brown Arguments of func: 3640873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3650873808bSJed Brown 3660873808bSJed Brown + adapt - adaptive controller context 3670873808bSJed Brown . ts - time stepping context 3680873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3690873808bSJed Brown 3700873808bSJed Brown Level: advanced 3710873808bSJed Brown 3720873808bSJed Brown .seealso: TSAdaptChoose() 3730873808bSJed Brown @*/ 374b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3750873808bSJed Brown { 3760873808bSJed Brown 3770873808bSJed Brown PetscFunctionBegin; 3780873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 37968ae941aSLisandro Dalcin adapt->checkstage = func; 3800873808bSJed Brown PetscFunctionReturn(0); 3810873808bSJed Brown } 3820873808bSJed Brown 3831c3436cfSJed Brown /*@ 384bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 385bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 386bf997491SLisandro Dalcin 3871917a363SLisandro Dalcin Logically collective on TSAdapt 388bf997491SLisandro Dalcin 389bf997491SLisandro Dalcin Input Arguments: 390bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 391bf997491SLisandro Dalcin - flag - whether to always accept steps 392bf997491SLisandro Dalcin 393bf997491SLisandro Dalcin Options Database Keys: 394ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 395bf997491SLisandro Dalcin 396bf997491SLisandro Dalcin Level: intermediate 397bf997491SLisandro Dalcin 398bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 399bf997491SLisandro Dalcin @*/ 400bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 401bf997491SLisandro Dalcin { 402bf997491SLisandro Dalcin PetscFunctionBegin; 403bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 404bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 405bf997491SLisandro Dalcin adapt->always_accept = flag; 406bf997491SLisandro Dalcin PetscFunctionReturn(0); 407bf997491SLisandro Dalcin } 408bf997491SLisandro Dalcin 409bf997491SLisandro Dalcin /*@ 4101917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4111c3436cfSJed Brown 4121917a363SLisandro Dalcin Logically collective on TSAdapt 4131917a363SLisandro Dalcin 4141917a363SLisandro Dalcin Input Arguments: 4151917a363SLisandro Dalcin + adapt - adaptive controller context 4161917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 417ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 4181917a363SLisandro Dalcin 4191917a363SLisandro Dalcin Options Database Keys: 420ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 421ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4221917a363SLisandro Dalcin 4231917a363SLisandro Dalcin Level: intermediate 4241917a363SLisandro Dalcin 4251917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4261917a363SLisandro Dalcin @*/ 4271917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4281917a363SLisandro Dalcin { 4291917a363SLisandro Dalcin PetscFunctionBegin; 4301917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4311917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4321917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4331917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4341917a363SLisandro 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); 4351917a363SLisandro 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); 4361917a363SLisandro 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); 4371917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4381917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4391917a363SLisandro Dalcin PetscFunctionReturn(0); 4401917a363SLisandro Dalcin } 4411917a363SLisandro Dalcin 4421917a363SLisandro Dalcin /*@ 4431917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4441917a363SLisandro Dalcin 4451917a363SLisandro Dalcin Not Collective 4461917a363SLisandro Dalcin 4471917a363SLisandro Dalcin Input Arguments: 4481917a363SLisandro Dalcin . adapt - adaptive controller context 4491917a363SLisandro Dalcin 4501917a363SLisandro Dalcin Ouput Arguments: 4511917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4521917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4531917a363SLisandro Dalcin 4541917a363SLisandro Dalcin Level: intermediate 4551917a363SLisandro Dalcin 4561917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4571917a363SLisandro Dalcin @*/ 4581917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4591917a363SLisandro Dalcin { 4601917a363SLisandro Dalcin PetscFunctionBegin; 4611917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4621917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4631917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4641917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4651917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4661917a363SLisandro Dalcin PetscFunctionReturn(0); 4671917a363SLisandro Dalcin } 4681917a363SLisandro Dalcin 4691917a363SLisandro Dalcin /*@ 4701917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4711917a363SLisandro Dalcin 4721917a363SLisandro Dalcin Logically collective on TSAdapt 4731917a363SLisandro Dalcin 4741917a363SLisandro Dalcin Input Arguments: 4751917a363SLisandro Dalcin + adapt - adaptive controller context 4761917a363SLisandro Dalcin . low - admissible decrease factor 477ec18b7bcSLisandro Dalcin - high - admissible increase factor 4781917a363SLisandro Dalcin 4791917a363SLisandro Dalcin Options Database Keys: 480ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 4811917a363SLisandro Dalcin 4821917a363SLisandro Dalcin Level: intermediate 4831917a363SLisandro Dalcin 4841917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptGetClip() 4851917a363SLisandro Dalcin @*/ 4861917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 4871917a363SLisandro Dalcin { 4881917a363SLisandro Dalcin PetscFunctionBegin; 4891917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4901917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 4911917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 4921917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 4931917a363SLisandro 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); 4941917a363SLisandro 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); 4951917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 4961917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 4971917a363SLisandro Dalcin PetscFunctionReturn(0); 4981917a363SLisandro Dalcin } 4991917a363SLisandro Dalcin 5001917a363SLisandro Dalcin /*@ 5011917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5021917a363SLisandro Dalcin 5031917a363SLisandro Dalcin Not Collective 5041917a363SLisandro Dalcin 5051917a363SLisandro Dalcin Input Arguments: 5061917a363SLisandro Dalcin . adapt - adaptive controller context 5071917a363SLisandro Dalcin 5081917a363SLisandro Dalcin Ouput Arguments: 5091917a363SLisandro Dalcin + low - optional, admissible decrease factor 5101917a363SLisandro Dalcin - high - optional, admissible increase factor 5111917a363SLisandro Dalcin 5121917a363SLisandro Dalcin Level: intermediate 5131917a363SLisandro Dalcin 5141917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptSetClip() 5151917a363SLisandro Dalcin @*/ 5161917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5171917a363SLisandro Dalcin { 5181917a363SLisandro Dalcin PetscFunctionBegin; 5191917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5201917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5211917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5221917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5231917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5241917a363SLisandro Dalcin PetscFunctionReturn(0); 5251917a363SLisandro Dalcin } 5261917a363SLisandro Dalcin 5271917a363SLisandro Dalcin /*@ 5281917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 5291917a363SLisandro Dalcin 5301917a363SLisandro Dalcin Logically collective on TSAdapt 5311c3436cfSJed Brown 5321c3436cfSJed Brown Input Arguments: 533552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5341c3436cfSJed Brown . hmin - minimum time step 5351c3436cfSJed Brown - hmax - maximum time step 5361c3436cfSJed Brown 5371c3436cfSJed Brown Options Database Keys: 538ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 539ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 5401c3436cfSJed Brown 5411c3436cfSJed Brown Level: intermediate 5421c3436cfSJed Brown 5431917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 5441c3436cfSJed Brown @*/ 5451c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 5461c3436cfSJed Brown { 5471c3436cfSJed Brown 5481c3436cfSJed Brown PetscFunctionBegin; 5494782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 550b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 551b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 552b1f5bccaSLisandro 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); 553b1f5bccaSLisandro 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); 554b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 555b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 556b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 557b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 558b1f5bccaSLisandro 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); 5591917a363SLisandro Dalcin PetscFunctionReturn(0); 5601917a363SLisandro Dalcin } 5611917a363SLisandro Dalcin 5621917a363SLisandro Dalcin /*@ 5631917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 5641917a363SLisandro Dalcin 5651917a363SLisandro Dalcin Not Collective 5661917a363SLisandro Dalcin 5671917a363SLisandro Dalcin Input Arguments: 5681917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5691917a363SLisandro Dalcin 5701917a363SLisandro Dalcin Output Arguments: 5711917a363SLisandro Dalcin + hmin - minimum time step 5721917a363SLisandro Dalcin - hmax - maximum time step 5731917a363SLisandro Dalcin 5741917a363SLisandro Dalcin Level: intermediate 5751917a363SLisandro Dalcin 5761917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 5771917a363SLisandro Dalcin @*/ 5781917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 5791917a363SLisandro Dalcin { 5801917a363SLisandro Dalcin 5811917a363SLisandro Dalcin PetscFunctionBegin; 5821917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5831917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 5841917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 5851917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 5861917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 5871c3436cfSJed Brown PetscFunctionReturn(0); 5881c3436cfSJed Brown } 5891c3436cfSJed Brown 590e55864a3SBarry Smith /* 59184df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 59284df9cb4SJed Brown 59384df9cb4SJed Brown Collective on TSAdapt 59484df9cb4SJed Brown 59584df9cb4SJed Brown Input Parameter: 59684df9cb4SJed Brown . adapt - the TSAdapt context 59784df9cb4SJed Brown 59884df9cb4SJed Brown Options Database Keys: 5991917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 600bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 6011917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6021917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6031917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 60423de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 60523de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 60623de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 607de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 608de50f1caSBarry 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 60984df9cb4SJed Brown 61084df9cb4SJed Brown Level: advanced 61184df9cb4SJed Brown 61284df9cb4SJed Brown Notes: 61384df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 61484df9cb4SJed Brown 61523de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 61684df9cb4SJed Brown 6171917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 6181917a363SLisandro Dalcin TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 619e55864a3SBarry Smith */ 6204416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 62184df9cb4SJed Brown { 62284df9cb4SJed Brown PetscErrorCode ierr; 62384df9cb4SJed Brown char type[256] = TSADAPTBASIC; 6241917a363SLisandro Dalcin PetscReal safety,reject_safety,clip[2],hmin,hmax; 6251c3436cfSJed Brown PetscBool set,flg; 6261917a363SLisandro Dalcin PetscInt two; 62784df9cb4SJed Brown 62884df9cb4SJed Brown PetscFunctionBegin; 6294782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 63084df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6311566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 632e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 63323de1d84SBarry 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); 63484df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 63584df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 63684df9cb4SJed Brown } 6371917a363SLisandro Dalcin 638bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 639bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 6401917a363SLisandro Dalcin 6411917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 6421917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 6431917a363SLisandro 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); 6441917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 6451917a363SLisandro Dalcin 6461917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 6471917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 6481917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 6491917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 6501917a363SLisandro Dalcin 6511917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 6521917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 6531917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 654bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 6551917a363SLisandro Dalcin 656*d580f011SEmil 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); 657*d580f011SEmil 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); 658*d580f011SEmil Constantinescu 6590298fd71SBarry 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); 6601917a363SLisandro Dalcin 6617619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 6627619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 6631917a363SLisandro Dalcin 664de50f1caSBarry 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); 665de50f1caSBarry Smith 6661917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 6671917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 6681917a363SLisandro Dalcin 669e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 67084df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 67184df9cb4SJed Brown PetscFunctionReturn(0); 67284df9cb4SJed Brown } 67384df9cb4SJed Brown 67484df9cb4SJed Brown /*@ 67584df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 67684df9cb4SJed Brown 6771917a363SLisandro Dalcin Logically collective on TSAdapt 67884df9cb4SJed Brown 67984df9cb4SJed Brown Input Argument: 68084df9cb4SJed Brown . adapt - adaptive controller 68184df9cb4SJed Brown 68284df9cb4SJed Brown Level: developer 68384df9cb4SJed Brown 68484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 68584df9cb4SJed Brown @*/ 68684df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 68784df9cb4SJed Brown { 68884df9cb4SJed Brown PetscErrorCode ierr; 68984df9cb4SJed Brown 69084df9cb4SJed Brown PetscFunctionBegin; 6914782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 69284df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 69384df9cb4SJed Brown PetscFunctionReturn(0); 69484df9cb4SJed Brown } 69584df9cb4SJed Brown 69684df9cb4SJed Brown /*@C 69784df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 69884df9cb4SJed Brown 6991917a363SLisandro Dalcin Logically collective on TSAdapt 70084df9cb4SJed Brown 70184df9cb4SJed Brown Input Arguments: 702552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 70384df9cb4SJed Brown . name - name of the candidate scheme to add 70484df9cb4SJed Brown . order - order of the candidate scheme 70584df9cb4SJed Brown . stageorder - stage order of the candidate scheme 7068d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 70784df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 70884df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 70984df9cb4SJed Brown 71084df9cb4SJed Brown Note: 71184df9cb4SJed Brown This routine is not available in Fortran. 71284df9cb4SJed Brown 71384df9cb4SJed Brown Level: developer 71484df9cb4SJed Brown 71584df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 71684df9cb4SJed Brown @*/ 7178d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 71884df9cb4SJed Brown { 71984df9cb4SJed Brown PetscInt c; 72084df9cb4SJed Brown 72184df9cb4SJed Brown PetscFunctionBegin; 72284df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 723ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 72484df9cb4SJed Brown if (inuse) { 725ce94432eSBarry 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()"); 72684df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 72784df9cb4SJed Brown } 7281c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7291c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 730bbd56ea5SKarl Rupp 73184df9cb4SJed Brown adapt->candidates.name[c] = name; 73284df9cb4SJed Brown adapt->candidates.order[c] = order; 73384df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7348d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 73584df9cb4SJed Brown adapt->candidates.cost[c] = cost; 73684df9cb4SJed Brown adapt->candidates.n++; 73784df9cb4SJed Brown PetscFunctionReturn(0); 73884df9cb4SJed Brown } 73984df9cb4SJed Brown 7408d59e960SJed Brown /*@C 7418d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 7428d59e960SJed Brown 7438d59e960SJed Brown Not Collective 7448d59e960SJed Brown 7458d59e960SJed Brown Input Arguments: 7468d59e960SJed Brown . adapt - time step adaptivity context 7478d59e960SJed Brown 7488d59e960SJed Brown Output Arguments: 7498d59e960SJed Brown + n - number of candidate schemes, always at least 1 7508d59e960SJed Brown . order - the order of each candidate scheme 7518d59e960SJed Brown . stageorder - the stage order of each candidate scheme 7528d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 7538d59e960SJed Brown - cost - the relative cost of each scheme 7548d59e960SJed Brown 7558d59e960SJed Brown Level: developer 7568d59e960SJed Brown 7578d59e960SJed Brown Note: 7588d59e960SJed Brown The current scheme is always returned in the first slot 7598d59e960SJed Brown 7608d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 7618d59e960SJed Brown @*/ 7628d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 7638d59e960SJed Brown { 7648d59e960SJed Brown PetscFunctionBegin; 7658d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7668d59e960SJed Brown if (n) *n = adapt->candidates.n; 7678d59e960SJed Brown if (order) *order = adapt->candidates.order; 7688d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 7698d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 7708d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 7718d59e960SJed Brown PetscFunctionReturn(0); 7728d59e960SJed Brown } 7738d59e960SJed Brown 77484df9cb4SJed Brown /*@C 77584df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 77684df9cb4SJed Brown 7771917a363SLisandro Dalcin Collective on TSAdapt 77884df9cb4SJed Brown 77984df9cb4SJed Brown Input Arguments: 78084df9cb4SJed Brown + adapt - adaptive contoller 78184df9cb4SJed Brown - h - current step size 78284df9cb4SJed Brown 78384df9cb4SJed Brown Output Arguments: 7841566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 78584df9cb4SJed Brown . next_h - step size to use for the next step 78684df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 78784df9cb4SJed Brown 7881c3436cfSJed Brown Note: 7891c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 7901c3436cfSJed Brown being retried after an initial rejection. 7911c3436cfSJed Brown 79284df9cb4SJed Brown Level: developer 79384df9cb4SJed Brown 79484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 79584df9cb4SJed Brown @*/ 79684df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 79784df9cb4SJed Brown { 79884df9cb4SJed Brown PetscErrorCode ierr; 7991566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 8001566a47fSLisandro Dalcin PetscInt scheme = 0; 8010b99f514SJed Brown PetscReal wlte = -1.0; 8025e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 8035e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 80484df9cb4SJed Brown 80584df9cb4SJed Brown PetscFunctionBegin; 80684df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 80784df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8081566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 80984df9cb4SJed Brown PetscValidPointer(next_h,5); 81084df9cb4SJed Brown PetscValidIntPointer(accept,6); 8111566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8121566a47fSLisandro Dalcin 8131566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8141566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8151566a47fSLisandro Dalcin *next_h = h; 8161566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8171566a47fSLisandro Dalcin PetscFunctionReturn(0); 8181566a47fSLisandro Dalcin } 8191566a47fSLisandro Dalcin 8205e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 8211566a47fSLisandro 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); 8221566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8231566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8241566a47fSLisandro Dalcin 82549354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 82636b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 82736b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 82836b54a69SLisandro Dalcin PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 829fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 830fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1]; 83136b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 832fe7350e0SStefano Zampini if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 83336b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 83449354f04SShri Abhyankar } 8351c3436cfSJed Brown 8361c3436cfSJed Brown if (adapt->monitor) { 8371566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8381c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8390b99f514SJed Brown if (wlte < 0) { 8400dea241dSBarry 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); 8410b99f514SJed Brown } else { 8420dea241dSBarry 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); 8430b99f514SJed Brown } 8441c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8451c3436cfSJed Brown } 84684df9cb4SJed Brown PetscFunctionReturn(0); 84784df9cb4SJed Brown } 84884df9cb4SJed Brown 84997335746SJed Brown /*@ 850de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 851de50f1caSBarry Smith before increasing the time step. 852de50f1caSBarry Smith 853de50f1caSBarry Smith Logicially Collective on TSAdapt 854de50f1caSBarry Smith 855de50f1caSBarry Smith Input Arguments: 856de50f1caSBarry Smith + adapt - adaptive controller context 857de50f1caSBarry Smith - cnt - the number of timesteps 858de50f1caSBarry Smith 859de50f1caSBarry Smith Options Database Key: 860de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 861de50f1caSBarry Smith 862de50f1caSBarry Smith Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 863de50f1caSBarry Smith The successful use of this option is problem dependent 864de50f1caSBarry Smith 865de50f1caSBarry Smith Developer Note: there is no theory to support this option 866de50f1caSBarry Smith 867de50f1caSBarry Smith Level: advanced 868de50f1caSBarry Smith 869de50f1caSBarry Smith .seealso: 870de50f1caSBarry Smith @*/ 871de50f1caSBarry Smith PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 872de50f1caSBarry Smith { 873de50f1caSBarry Smith PetscFunctionBegin; 874de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt; 875de50f1caSBarry Smith PetscFunctionReturn(0); 876de50f1caSBarry Smith } 877de50f1caSBarry Smith 878de50f1caSBarry Smith 879de50f1caSBarry Smith /*@ 88097335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 88197335746SJed Brown 8821917a363SLisandro Dalcin Collective on TSAdapt 88397335746SJed Brown 88497335746SJed Brown Input Arguments: 88597335746SJed Brown + adapt - adaptive controller context 886b295832fSPierre Barbier de Reuille . ts - time stepper 887b295832fSPierre Barbier de Reuille . t - Current simulation time 888b295832fSPierre Barbier de Reuille - Y - Current solution vector 88997335746SJed Brown 89097335746SJed Brown Output Arguments: 89197335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 89297335746SJed Brown 89397335746SJed Brown Level: developer 89497335746SJed Brown 89597335746SJed Brown .seealso: 89697335746SJed Brown @*/ 897b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 89897335746SJed Brown { 89997335746SJed Brown PetscErrorCode ierr; 9001566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 90197335746SJed Brown 90297335746SJed Brown PetscFunctionBegin; 9034782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 9044782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 9054782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 9061566a47fSLisandro Dalcin 9071566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 90897335746SJed Brown if (snesreason < 0) { 90997335746SJed Brown *accept = PETSC_FALSE; 9106de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 91197335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 91248c19aefSLisandro 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); 91397335746SJed Brown if (adapt->monitor) { 91497335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 9150dea241dSBarry 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); 91697335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 91797335746SJed Brown } 918cb9d8021SPierre Barbier de Reuille } 919cb9d8021SPierre Barbier de Reuille } else { 9201566a47fSLisandro Dalcin *accept = PETSC_TRUE; 921b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 922cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 923b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 924cb9d8021SPierre Barbier de Reuille } 925cb9d8021SPierre Barbier de Reuille } 926cb9d8021SPierre Barbier de Reuille 9271566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 9281566a47fSLisandro Dalcin PetscReal dt,new_dt; 9291566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 930cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 93197335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 932de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 93397335746SJed Brown if (adapt->monitor) { 93497335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 935e6d0a238SBarry 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); 93697335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 93797335746SJed Brown } 93897335746SJed Brown } 93997335746SJed Brown PetscFunctionReturn(0); 94097335746SJed Brown } 94197335746SJed Brown 94284df9cb4SJed Brown /*@ 94384df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 94484df9cb4SJed Brown 94584df9cb4SJed Brown Collective on MPI_Comm 94684df9cb4SJed Brown 94784df9cb4SJed Brown Input Parameter: 94884df9cb4SJed Brown . comm - The communicator 94984df9cb4SJed Brown 95084df9cb4SJed Brown Output Parameter: 95184df9cb4SJed Brown . adapt - new TSAdapt object 95284df9cb4SJed Brown 95384df9cb4SJed Brown Level: developer 95484df9cb4SJed Brown 95584df9cb4SJed Brown Notes: 95684df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 95784df9cb4SJed Brown 95884df9cb4SJed Brown .keywords: TSAdapt, create 959552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 96084df9cb4SJed Brown @*/ 96184df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 96284df9cb4SJed Brown { 96384df9cb4SJed Brown PetscErrorCode ierr; 96484df9cb4SJed Brown TSAdapt adapt; 96584df9cb4SJed Brown 96684df9cb4SJed Brown PetscFunctionBegin; 9673b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 9683b3bcf4cSLisandro Dalcin *inadapt = NULL; 9693b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 9703b3bcf4cSLisandro Dalcin 97173107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 9721c3436cfSJed Brown 973bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 9741917a363SLisandro Dalcin adapt->safety = 0.9; 9751917a363SLisandro Dalcin adapt->reject_safety = 0.5; 9761917a363SLisandro Dalcin adapt->clip[0] = 0.1; 9771917a363SLisandro Dalcin adapt->clip[1] = 10.; 9781c3436cfSJed Brown adapt->dt_min = 1e-20; 9791917a363SLisandro Dalcin adapt->dt_max = 1e+20; 980*d580f011SEmil Constantinescu adapt->ignore_max = -PETSC_INFINITY; 981*d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE; 98297335746SJed Brown adapt->scale_solve_failed = 0.25; 983fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 984fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 985fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 986fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 9877619abb3SShri adapt->wnormtype = NORM_2; 988de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0; 9891917a363SLisandro Dalcin 99084df9cb4SJed Brown *inadapt = adapt; 99184df9cb4SJed Brown PetscFunctionReturn(0); 99284df9cb4SJed Brown } 993