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); 1584df9cb4SJed Brown 1684df9cb4SJed Brown /*@C 171c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation 181c84c290SBarry Smith 191c84c290SBarry Smith Not Collective 201c84c290SBarry Smith 211c84c290SBarry Smith Input Parameters: 221c84c290SBarry Smith + name_scheme - name of user-defined adaptivity scheme 231c84c290SBarry Smith - routine_create - routine to create method context 241c84c290SBarry Smith 251c84c290SBarry Smith Notes: 261c84c290SBarry Smith TSAdaptRegister() may be called multiple times to add several user-defined families. 271c84c290SBarry Smith 281c84c290SBarry Smith Sample usage: 291c84c290SBarry Smith .vb 30bdf89e91SBarry Smith TSAdaptRegister("my_scheme",MySchemeCreate); 311c84c290SBarry Smith .ve 321c84c290SBarry Smith 331c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via 341c84c290SBarry Smith $ TSAdaptSetType(ts,"my_scheme") 351c84c290SBarry Smith or at runtime via the option 361c84c290SBarry Smith $ -ts_adapt_type my_scheme 3784df9cb4SJed Brown 3884df9cb4SJed Brown Level: advanced 391c84c290SBarry Smith 401c84c290SBarry Smith .keywords: TSAdapt, register 411c84c290SBarry Smith 421c84c290SBarry Smith .seealso: TSAdaptRegisterAll() 4384df9cb4SJed Brown @*/ 44bdf89e91SBarry Smith PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 4584df9cb4SJed Brown { 4684df9cb4SJed Brown PetscErrorCode ierr; 4784df9cb4SJed Brown 4884df9cb4SJed Brown PetscFunctionBegin; 49a240a19fSJed Brown ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 5084df9cb4SJed Brown PetscFunctionReturn(0); 5184df9cb4SJed Brown } 5284df9cb4SJed Brown 5384df9cb4SJed Brown /*@C 5484df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5584df9cb4SJed Brown 5684df9cb4SJed Brown Not Collective 5784df9cb4SJed Brown 5884df9cb4SJed Brown Level: advanced 5984df9cb4SJed Brown 6084df9cb4SJed Brown .keywords: TSAdapt, register, all 6184df9cb4SJed Brown 6284df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6384df9cb4SJed Brown @*/ 64607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6584df9cb4SJed Brown { 6684df9cb4SJed Brown PetscErrorCode ierr; 6784df9cb4SJed Brown 6884df9cb4SJed Brown PetscFunctionBegin; 690f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 700f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 71bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 721566a47fSLisandro Dalcin ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 73cb7ab186SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP);CHKERRQ(ierr); 74bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 751917a363SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE);CHKERRQ(ierr); 7684df9cb4SJed Brown PetscFunctionReturn(0); 7784df9cb4SJed Brown } 7884df9cb4SJed Brown 7984df9cb4SJed Brown /*@C 80560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 8184df9cb4SJed Brown called from PetscFinalize(). 8284df9cb4SJed Brown 8384df9cb4SJed Brown Level: developer 8484df9cb4SJed Brown 8584df9cb4SJed Brown .keywords: Petsc, destroy, package 8684df9cb4SJed Brown .seealso: PetscFinalize() 8784df9cb4SJed Brown @*/ 8884df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8984df9cb4SJed Brown { 9037e93019SBarry Smith PetscErrorCode ierr; 9137e93019SBarry Smith 9284df9cb4SJed Brown PetscFunctionBegin; 9337e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9484df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9584df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9684df9cb4SJed Brown PetscFunctionReturn(0); 9784df9cb4SJed Brown } 9884df9cb4SJed Brown 9984df9cb4SJed Brown /*@C 10084df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 10184df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 1021917a363SLisandro Dalcin TSAdaptCreate() when using static libraries. 10384df9cb4SJed Brown 10484df9cb4SJed Brown Level: developer 10584df9cb4SJed Brown 10684df9cb4SJed Brown .keywords: TSAdapt, initialize, package 10784df9cb4SJed Brown .seealso: PetscInitialize() 10884df9cb4SJed Brown @*/ 109607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11084df9cb4SJed Brown { 11184df9cb4SJed Brown PetscErrorCode ierr; 11284df9cb4SJed Brown 11384df9cb4SJed Brown PetscFunctionBegin; 11484df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11584df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11684df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 117607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 11884df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 11984df9cb4SJed Brown PetscFunctionReturn(0); 12084df9cb4SJed Brown } 12184df9cb4SJed Brown 1227eef6055SBarry Smith /*@C 1237eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1247eef6055SBarry Smith 1257eef6055SBarry Smith Logicially Collective on TSAdapt 1267eef6055SBarry Smith 1277eef6055SBarry Smith Input Parameter: 128d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1297eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1307eef6055SBarry Smith 1317eef6055SBarry Smith Options Database: 132*ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type 1337eef6055SBarry Smith 1347eef6055SBarry Smith Level: intermediate 1357eef6055SBarry Smith 1367eef6055SBarry Smith .keywords: TSAdapt, create 1377eef6055SBarry Smith 1387eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1397eef6055SBarry Smith @*/ 14019fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 14184df9cb4SJed Brown { 142ef749922SLisandro Dalcin PetscBool match; 14384df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 14484df9cb4SJed Brown 14584df9cb4SJed Brown PetscFunctionBegin; 1464782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 147b92453a8SLisandro Dalcin PetscValidCharPointer(type,2); 148ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 149ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1501c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 15184df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 152ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1534cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 15484df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15568ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15684df9cb4SJed Brown PetscFunctionReturn(0); 15784df9cb4SJed Brown } 15884df9cb4SJed Brown 159d0288e4fSLisandro Dalcin /*@C 160d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 161d0288e4fSLisandro Dalcin 162d0288e4fSLisandro Dalcin Not Collective 163d0288e4fSLisandro Dalcin 164d0288e4fSLisandro Dalcin Input Parameter: 165d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 166d0288e4fSLisandro Dalcin 167d0288e4fSLisandro Dalcin Output Parameter: 168d0288e4fSLisandro Dalcin . type - The name of TS adapter method 169d0288e4fSLisandro Dalcin 170d0288e4fSLisandro Dalcin Level: intermediate 171d0288e4fSLisandro Dalcin 172d0288e4fSLisandro Dalcin .keywords: TSAdapt, get, type 173d0288e4fSLisandro Dalcin .seealso TSAdaptSetType() 174d0288e4fSLisandro Dalcin @*/ 175d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 176d0288e4fSLisandro Dalcin { 177d0288e4fSLisandro Dalcin PetscFunctionBegin; 178d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 179d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 180d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 181d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 182d0288e4fSLisandro Dalcin } 183d0288e4fSLisandro Dalcin 18484df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 18584df9cb4SJed Brown { 18684df9cb4SJed Brown PetscErrorCode ierr; 18784df9cb4SJed Brown 18884df9cb4SJed Brown PetscFunctionBegin; 1894782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 19084df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 19184df9cb4SJed Brown PetscFunctionReturn(0); 19284df9cb4SJed Brown } 19384df9cb4SJed Brown 194ad6bc421SBarry Smith /*@C 195ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 196ad6bc421SBarry Smith 197ad6bc421SBarry Smith Collective on PetscViewer 198ad6bc421SBarry Smith 199ad6bc421SBarry Smith Input Parameters: 200ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 201ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 202ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 203ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 204ad6bc421SBarry Smith 205ad6bc421SBarry Smith Level: intermediate 206ad6bc421SBarry Smith 207ad6bc421SBarry Smith Notes: 208ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 209ad6bc421SBarry Smith 210ad6bc421SBarry Smith Notes for advanced users: 211ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 212ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 213ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 214ad6bc421SBarry Smith format is 215ad6bc421SBarry Smith .vb 216ad6bc421SBarry Smith has not yet been determined 217ad6bc421SBarry Smith .ve 218ad6bc421SBarry Smith 219ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 220ad6bc421SBarry Smith @*/ 2214782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 222ad6bc421SBarry Smith { 223ad6bc421SBarry Smith PetscErrorCode ierr; 224ad6bc421SBarry Smith PetscBool isbinary; 225ad6bc421SBarry Smith char type[256]; 226ad6bc421SBarry Smith 227ad6bc421SBarry Smith PetscFunctionBegin; 2284782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 229ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 230ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 231ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 232ad6bc421SBarry Smith 233060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2344782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2354782b174SLisandro Dalcin if (adapt->ops->load) { 2364782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 237ad6bc421SBarry Smith } 238ad6bc421SBarry Smith PetscFunctionReturn(0); 239ad6bc421SBarry Smith } 240ad6bc421SBarry Smith 24184df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 24284df9cb4SJed Brown { 24384df9cb4SJed Brown PetscErrorCode ierr; 2441917a363SLisandro Dalcin PetscBool iascii,isbinary,isnone; 24584df9cb4SJed Brown 24684df9cb4SJed Brown PetscFunctionBegin; 2474782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2484782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2494782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2504782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 251251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 252ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 25384df9cb4SJed Brown if (iascii) { 254dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 2551917a363SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr); 2561917a363SLisandro Dalcin if (!isnone) { 257bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 2581917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 2591917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 2601917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 2611917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 262bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 263bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 2641917a363SLisandro Dalcin } 26584df9cb4SJed Brown if (adapt->ops->view) { 26684df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 26784df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 26884df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 26984df9cb4SJed Brown } 270ad6bc421SBarry Smith } else if (isbinary) { 271ad6bc421SBarry Smith char type[256]; 272ad6bc421SBarry Smith 273ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 274ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 275ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 276bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 277f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 278f2c2a1b9SBarry Smith } 27984df9cb4SJed Brown PetscFunctionReturn(0); 28084df9cb4SJed Brown } 28184df9cb4SJed Brown 282099af0a3SLisandro Dalcin /*@ 283099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 284099af0a3SLisandro Dalcin 285099af0a3SLisandro Dalcin Collective on TS 286099af0a3SLisandro Dalcin 287099af0a3SLisandro Dalcin Input Parameter: 288099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 289099af0a3SLisandro Dalcin 290099af0a3SLisandro Dalcin Level: developer 291099af0a3SLisandro Dalcin 292099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 293099af0a3SLisandro Dalcin @*/ 294099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 295099af0a3SLisandro Dalcin { 296099af0a3SLisandro Dalcin PetscErrorCode ierr; 297099af0a3SLisandro Dalcin 298099af0a3SLisandro Dalcin PetscFunctionBegin; 299099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 300099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 301099af0a3SLisandro Dalcin PetscFunctionReturn(0); 302099af0a3SLisandro Dalcin } 303099af0a3SLisandro Dalcin 30484df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 30584df9cb4SJed Brown { 30684df9cb4SJed Brown PetscErrorCode ierr; 30784df9cb4SJed Brown 30884df9cb4SJed Brown PetscFunctionBegin; 30984df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 31084df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3114782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 312099af0a3SLisandro Dalcin 313099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 314099af0a3SLisandro Dalcin 31584df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3161c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 31784df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 31884df9cb4SJed Brown PetscFunctionReturn(0); 31984df9cb4SJed Brown } 32084df9cb4SJed Brown 3211c3436cfSJed Brown /*@ 3221c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3231c3436cfSJed Brown 3241c3436cfSJed Brown Collective on TSAdapt 3251c3436cfSJed Brown 3261c3436cfSJed Brown Input Arguments: 3271c3436cfSJed Brown + adapt - adaptive controller context 3281c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3291c3436cfSJed Brown 330bf997491SLisandro Dalcin Options Database Keys: 331*ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring 332bf997491SLisandro Dalcin 3331c3436cfSJed Brown Level: intermediate 3341c3436cfSJed Brown 3351c3436cfSJed Brown .seealso: TSAdaptChoose() 3361c3436cfSJed Brown @*/ 3371c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3381c3436cfSJed Brown { 3391c3436cfSJed Brown PetscErrorCode ierr; 3401c3436cfSJed Brown 3411c3436cfSJed Brown PetscFunctionBegin; 3424782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3434782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3441c3436cfSJed Brown if (flg) { 345ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3461c3436cfSJed Brown } else { 3471c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3481c3436cfSJed Brown } 3491c3436cfSJed Brown PetscFunctionReturn(0); 3501c3436cfSJed Brown } 3511c3436cfSJed Brown 3520873808bSJed Brown /*@C 353bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3540873808bSJed Brown 3550873808bSJed Brown Logically collective on TSAdapt 3560873808bSJed Brown 3570873808bSJed Brown Input Arguments: 3580873808bSJed Brown + adapt - adaptive controller context 3590873808bSJed Brown - func - stage check function 3600873808bSJed Brown 3610873808bSJed Brown Arguments of func: 3620873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3630873808bSJed Brown 3640873808bSJed Brown + adapt - adaptive controller context 3650873808bSJed Brown . ts - time stepping context 3660873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3670873808bSJed Brown 3680873808bSJed Brown Level: advanced 3690873808bSJed Brown 3700873808bSJed Brown .seealso: TSAdaptChoose() 3710873808bSJed Brown @*/ 372b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3730873808bSJed Brown { 3740873808bSJed Brown 3750873808bSJed Brown PetscFunctionBegin; 3760873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 37768ae941aSLisandro Dalcin adapt->checkstage = func; 3780873808bSJed Brown PetscFunctionReturn(0); 3790873808bSJed Brown } 3800873808bSJed Brown 3811c3436cfSJed Brown /*@ 382bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 383bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 384bf997491SLisandro Dalcin 3851917a363SLisandro Dalcin Logically collective on TSAdapt 386bf997491SLisandro Dalcin 387bf997491SLisandro Dalcin Input Arguments: 388bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 389bf997491SLisandro Dalcin - flag - whether to always accept steps 390bf997491SLisandro Dalcin 391bf997491SLisandro Dalcin Options Database Keys: 392*ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps 393bf997491SLisandro Dalcin 394bf997491SLisandro Dalcin Level: intermediate 395bf997491SLisandro Dalcin 396bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 397bf997491SLisandro Dalcin @*/ 398bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 399bf997491SLisandro Dalcin { 400bf997491SLisandro Dalcin PetscFunctionBegin; 401bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 402bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 403bf997491SLisandro Dalcin adapt->always_accept = flag; 404bf997491SLisandro Dalcin PetscFunctionReturn(0); 405bf997491SLisandro Dalcin } 406bf997491SLisandro Dalcin 407bf997491SLisandro Dalcin /*@ 4081917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4091c3436cfSJed Brown 4101917a363SLisandro Dalcin Logically collective on TSAdapt 4111917a363SLisandro Dalcin 4121917a363SLisandro Dalcin Input Arguments: 4131917a363SLisandro Dalcin + adapt - adaptive controller context 4141917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 415*ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected 4161917a363SLisandro Dalcin 4171917a363SLisandro Dalcin Options Database Keys: 418*ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor 419*ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 4201917a363SLisandro Dalcin 4211917a363SLisandro Dalcin Level: intermediate 4221917a363SLisandro Dalcin 4231917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 4241917a363SLisandro Dalcin @*/ 4251917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 4261917a363SLisandro Dalcin { 4271917a363SLisandro Dalcin PetscFunctionBegin; 4281917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4291917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 4301917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 4311917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 4321917a363SLisandro 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); 4331917a363SLisandro 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); 4341917a363SLisandro 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); 4351917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 4361917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 4371917a363SLisandro Dalcin PetscFunctionReturn(0); 4381917a363SLisandro Dalcin } 4391917a363SLisandro Dalcin 4401917a363SLisandro Dalcin /*@ 4411917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 4421917a363SLisandro Dalcin 4431917a363SLisandro Dalcin Not Collective 4441917a363SLisandro Dalcin 4451917a363SLisandro Dalcin Input Arguments: 4461917a363SLisandro Dalcin . adapt - adaptive controller context 4471917a363SLisandro Dalcin 4481917a363SLisandro Dalcin Ouput Arguments: 4491917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 4501917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 4511917a363SLisandro Dalcin 4521917a363SLisandro Dalcin Level: intermediate 4531917a363SLisandro Dalcin 4541917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 4551917a363SLisandro Dalcin @*/ 4561917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 4571917a363SLisandro Dalcin { 4581917a363SLisandro Dalcin PetscFunctionBegin; 4591917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4601917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 4611917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 4621917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 4631917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 4641917a363SLisandro Dalcin PetscFunctionReturn(0); 4651917a363SLisandro Dalcin } 4661917a363SLisandro Dalcin 4671917a363SLisandro Dalcin /*@ 4681917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 4691917a363SLisandro Dalcin 4701917a363SLisandro Dalcin Logically collective on TSAdapt 4711917a363SLisandro Dalcin 4721917a363SLisandro Dalcin Input Arguments: 4731917a363SLisandro Dalcin + adapt - adaptive controller context 4741917a363SLisandro Dalcin . low - admissible decrease factor 475*ec18b7bcSLisandro Dalcin - high - admissible increase factor 4761917a363SLisandro Dalcin 4771917a363SLisandro Dalcin Options Database Keys: 478*ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 4791917a363SLisandro Dalcin 4801917a363SLisandro Dalcin Level: intermediate 4811917a363SLisandro Dalcin 4821917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptGetClip() 4831917a363SLisandro Dalcin @*/ 4841917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 4851917a363SLisandro Dalcin { 4861917a363SLisandro Dalcin PetscFunctionBegin; 4871917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4881917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 4891917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 4901917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 4911917a363SLisandro 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); 4921917a363SLisandro 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); 4931917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 4941917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 4951917a363SLisandro Dalcin PetscFunctionReturn(0); 4961917a363SLisandro Dalcin } 4971917a363SLisandro Dalcin 4981917a363SLisandro Dalcin /*@ 4991917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 5001917a363SLisandro Dalcin 5011917a363SLisandro Dalcin Not Collective 5021917a363SLisandro Dalcin 5031917a363SLisandro Dalcin Input Arguments: 5041917a363SLisandro Dalcin . adapt - adaptive controller context 5051917a363SLisandro Dalcin 5061917a363SLisandro Dalcin Ouput Arguments: 5071917a363SLisandro Dalcin + low - optional, admissible decrease factor 5081917a363SLisandro Dalcin - high - optional, admissible increase factor 5091917a363SLisandro Dalcin 5101917a363SLisandro Dalcin Level: intermediate 5111917a363SLisandro Dalcin 5121917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptSetClip() 5131917a363SLisandro Dalcin @*/ 5141917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 5151917a363SLisandro Dalcin { 5161917a363SLisandro Dalcin PetscFunctionBegin; 5171917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5181917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 5191917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 5201917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 5211917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 5221917a363SLisandro Dalcin PetscFunctionReturn(0); 5231917a363SLisandro Dalcin } 5241917a363SLisandro Dalcin 5251917a363SLisandro Dalcin /*@ 5261917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 5271917a363SLisandro Dalcin 5281917a363SLisandro Dalcin Logically collective on TSAdapt 5291c3436cfSJed Brown 5301c3436cfSJed Brown Input Arguments: 531552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5321c3436cfSJed Brown . hmin - minimum time step 5331c3436cfSJed Brown - hmax - maximum time step 5341c3436cfSJed Brown 5351c3436cfSJed Brown Options Database Keys: 536*ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step 537*ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step 5381c3436cfSJed Brown 5391c3436cfSJed Brown Level: intermediate 5401c3436cfSJed Brown 5411917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 5421c3436cfSJed Brown @*/ 5431c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 5441c3436cfSJed Brown { 5451c3436cfSJed Brown 5461c3436cfSJed Brown PetscFunctionBegin; 5474782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 548b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 549b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 550b1f5bccaSLisandro 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); 551b1f5bccaSLisandro 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); 552b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 553b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 554b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 555b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 556b1f5bccaSLisandro 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); 5571917a363SLisandro Dalcin PetscFunctionReturn(0); 5581917a363SLisandro Dalcin } 5591917a363SLisandro Dalcin 5601917a363SLisandro Dalcin /*@ 5611917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 5621917a363SLisandro Dalcin 5631917a363SLisandro Dalcin Not Collective 5641917a363SLisandro Dalcin 5651917a363SLisandro Dalcin Input Arguments: 5661917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5671917a363SLisandro Dalcin 5681917a363SLisandro Dalcin Output Arguments: 5691917a363SLisandro Dalcin + hmin - minimum time step 5701917a363SLisandro Dalcin - hmax - maximum time step 5711917a363SLisandro Dalcin 5721917a363SLisandro Dalcin Level: intermediate 5731917a363SLisandro Dalcin 5741917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 5751917a363SLisandro Dalcin @*/ 5761917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 5771917a363SLisandro Dalcin { 5781917a363SLisandro Dalcin 5791917a363SLisandro Dalcin PetscFunctionBegin; 5801917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5811917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 5821917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 5831917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 5841917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 5851c3436cfSJed Brown PetscFunctionReturn(0); 5861c3436cfSJed Brown } 5871c3436cfSJed Brown 588e55864a3SBarry Smith /* 58984df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 59084df9cb4SJed Brown 59184df9cb4SJed Brown Collective on TSAdapt 59284df9cb4SJed Brown 59384df9cb4SJed Brown Input Parameter: 59484df9cb4SJed Brown . adapt - the TSAdapt context 59584df9cb4SJed Brown 59684df9cb4SJed Brown Options Database Keys: 5971917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 598bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 5991917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 6001917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 6011917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 60223de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 60323de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 60423de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 60523de1d84SBarry Smith - -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 60684df9cb4SJed Brown 60784df9cb4SJed Brown Level: advanced 60884df9cb4SJed Brown 60984df9cb4SJed Brown Notes: 61084df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 61184df9cb4SJed Brown 61223de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 61384df9cb4SJed Brown 6141917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 6151917a363SLisandro Dalcin TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 616e55864a3SBarry Smith */ 6174416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 61884df9cb4SJed Brown { 61984df9cb4SJed Brown PetscErrorCode ierr; 62084df9cb4SJed Brown char type[256] = TSADAPTBASIC; 6211917a363SLisandro Dalcin PetscReal safety,reject_safety,clip[2],hmin,hmax; 6221c3436cfSJed Brown PetscBool set,flg; 6231917a363SLisandro Dalcin PetscInt two; 62484df9cb4SJed Brown 62584df9cb4SJed Brown PetscFunctionBegin; 6264782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 62784df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6281566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 629e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 63023de1d84SBarry 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); 63184df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 63284df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 63384df9cb4SJed Brown } 6341917a363SLisandro Dalcin 635bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 636bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 6371917a363SLisandro Dalcin 6381917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 6391917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 6401917a363SLisandro 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); 6411917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 6421917a363SLisandro Dalcin 6431917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 6441917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 6451917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 6461917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 6471917a363SLisandro Dalcin 6481917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 6491917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 6501917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 651bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 6521917a363SLisandro Dalcin 6530298fd71SBarry 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); 6541917a363SLisandro Dalcin 6557619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 6567619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 6571917a363SLisandro Dalcin 6581917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 6591917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 6601917a363SLisandro Dalcin 661e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 66284df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 66384df9cb4SJed Brown PetscFunctionReturn(0); 66484df9cb4SJed Brown } 66584df9cb4SJed Brown 66684df9cb4SJed Brown /*@ 66784df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 66884df9cb4SJed Brown 6691917a363SLisandro Dalcin Logically collective on TSAdapt 67084df9cb4SJed Brown 67184df9cb4SJed Brown Input Argument: 67284df9cb4SJed Brown . adapt - adaptive controller 67384df9cb4SJed Brown 67484df9cb4SJed Brown Level: developer 67584df9cb4SJed Brown 67684df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 67784df9cb4SJed Brown @*/ 67884df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 67984df9cb4SJed Brown { 68084df9cb4SJed Brown PetscErrorCode ierr; 68184df9cb4SJed Brown 68284df9cb4SJed Brown PetscFunctionBegin; 6834782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 68484df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 68584df9cb4SJed Brown PetscFunctionReturn(0); 68684df9cb4SJed Brown } 68784df9cb4SJed Brown 68884df9cb4SJed Brown /*@C 68984df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 69084df9cb4SJed Brown 6911917a363SLisandro Dalcin Logically collective on TSAdapt 69284df9cb4SJed Brown 69384df9cb4SJed Brown Input Arguments: 694552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 69584df9cb4SJed Brown . name - name of the candidate scheme to add 69684df9cb4SJed Brown . order - order of the candidate scheme 69784df9cb4SJed Brown . stageorder - stage order of the candidate scheme 6988d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 69984df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 70084df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 70184df9cb4SJed Brown 70284df9cb4SJed Brown Note: 70384df9cb4SJed Brown This routine is not available in Fortran. 70484df9cb4SJed Brown 70584df9cb4SJed Brown Level: developer 70684df9cb4SJed Brown 70784df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 70884df9cb4SJed Brown @*/ 7098d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 71084df9cb4SJed Brown { 71184df9cb4SJed Brown PetscInt c; 71284df9cb4SJed Brown 71384df9cb4SJed Brown PetscFunctionBegin; 71484df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 715ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 71684df9cb4SJed Brown if (inuse) { 717ce94432eSBarry 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()"); 71884df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 71984df9cb4SJed Brown } 7201c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7211c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 722bbd56ea5SKarl Rupp 72384df9cb4SJed Brown adapt->candidates.name[c] = name; 72484df9cb4SJed Brown adapt->candidates.order[c] = order; 72584df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7268d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 72784df9cb4SJed Brown adapt->candidates.cost[c] = cost; 72884df9cb4SJed Brown adapt->candidates.n++; 72984df9cb4SJed Brown PetscFunctionReturn(0); 73084df9cb4SJed Brown } 73184df9cb4SJed Brown 7328d59e960SJed Brown /*@C 7338d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 7348d59e960SJed Brown 7358d59e960SJed Brown Not Collective 7368d59e960SJed Brown 7378d59e960SJed Brown Input Arguments: 7388d59e960SJed Brown . adapt - time step adaptivity context 7398d59e960SJed Brown 7408d59e960SJed Brown Output Arguments: 7418d59e960SJed Brown + n - number of candidate schemes, always at least 1 7428d59e960SJed Brown . order - the order of each candidate scheme 7438d59e960SJed Brown . stageorder - the stage order of each candidate scheme 7448d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 7458d59e960SJed Brown - cost - the relative cost of each scheme 7468d59e960SJed Brown 7478d59e960SJed Brown Level: developer 7488d59e960SJed Brown 7498d59e960SJed Brown Note: 7508d59e960SJed Brown The current scheme is always returned in the first slot 7518d59e960SJed Brown 7528d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 7538d59e960SJed Brown @*/ 7548d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 7558d59e960SJed Brown { 7568d59e960SJed Brown PetscFunctionBegin; 7578d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7588d59e960SJed Brown if (n) *n = adapt->candidates.n; 7598d59e960SJed Brown if (order) *order = adapt->candidates.order; 7608d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 7618d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 7628d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 7638d59e960SJed Brown PetscFunctionReturn(0); 7648d59e960SJed Brown } 7658d59e960SJed Brown 76684df9cb4SJed Brown /*@C 76784df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 76884df9cb4SJed Brown 7691917a363SLisandro Dalcin Collective on TSAdapt 77084df9cb4SJed Brown 77184df9cb4SJed Brown Input Arguments: 77284df9cb4SJed Brown + adapt - adaptive contoller 77384df9cb4SJed Brown - h - current step size 77484df9cb4SJed Brown 77584df9cb4SJed Brown Output Arguments: 7761566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 77784df9cb4SJed Brown . next_h - step size to use for the next step 77884df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 77984df9cb4SJed Brown 7801c3436cfSJed Brown Note: 7811c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 7821c3436cfSJed Brown being retried after an initial rejection. 7831c3436cfSJed Brown 78484df9cb4SJed Brown Level: developer 78584df9cb4SJed Brown 78684df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 78784df9cb4SJed Brown @*/ 78884df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 78984df9cb4SJed Brown { 79084df9cb4SJed Brown PetscErrorCode ierr; 7911566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 7921566a47fSLisandro Dalcin PetscInt scheme = 0; 7930b99f514SJed Brown PetscReal wlte = -1.0; 7945e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 7955e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 79684df9cb4SJed Brown 79784df9cb4SJed Brown PetscFunctionBegin; 79884df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 79984df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8001566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 80184df9cb4SJed Brown PetscValidPointer(next_h,5); 80284df9cb4SJed Brown PetscValidIntPointer(accept,6); 8031566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8041566a47fSLisandro Dalcin 8051566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8061566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8071566a47fSLisandro Dalcin *next_h = h; 8081566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8091566a47fSLisandro Dalcin PetscFunctionReturn(0); 8101566a47fSLisandro Dalcin } 8111566a47fSLisandro Dalcin 8125e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 8131566a47fSLisandro 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); 8141566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8151566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8161566a47fSLisandro Dalcin 81749354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 81836b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 81936b54a69SLisandro Dalcin PetscReal t = ts->ptime + ts->time_step, h = *next_h; 82036b54a69SLisandro Dalcin PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 82136b54a69SLisandro Dalcin PetscReal a = (PetscReal)1.01; /* allow 1% step size increase */ 82236b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax; 82336b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h > hmax/2) *next_h = hmax/2; 82436b54a69SLisandro Dalcin if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 82549354f04SShri Abhyankar } 8261c3436cfSJed Brown 8271c3436cfSJed Brown if (adapt->monitor) { 8281566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8291c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8300b99f514SJed Brown if (wlte < 0) { 8310dea241dSBarry 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); 8320b99f514SJed Brown } else { 8330dea241dSBarry 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); 8340b99f514SJed Brown } 8351c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8361c3436cfSJed Brown } 83784df9cb4SJed Brown PetscFunctionReturn(0); 83884df9cb4SJed Brown } 83984df9cb4SJed Brown 84097335746SJed Brown /*@ 84197335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 84297335746SJed Brown 8431917a363SLisandro Dalcin Collective on TSAdapt 84497335746SJed Brown 84597335746SJed Brown Input Arguments: 84697335746SJed Brown + adapt - adaptive controller context 847b295832fSPierre Barbier de Reuille . ts - time stepper 848b295832fSPierre Barbier de Reuille . t - Current simulation time 849b295832fSPierre Barbier de Reuille - Y - Current solution vector 85097335746SJed Brown 85197335746SJed Brown Output Arguments: 85297335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 85397335746SJed Brown 85497335746SJed Brown Level: developer 85597335746SJed Brown 85697335746SJed Brown .seealso: 85797335746SJed Brown @*/ 858b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 85997335746SJed Brown { 86097335746SJed Brown PetscErrorCode ierr; 8611566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 86297335746SJed Brown 86397335746SJed Brown PetscFunctionBegin; 8644782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8654782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8664782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 8671566a47fSLisandro Dalcin 8681566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 86997335746SJed Brown if (snesreason < 0) { 87097335746SJed Brown *accept = PETSC_FALSE; 8716de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 87297335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 87348c19aefSLisandro 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); 87497335746SJed Brown if (adapt->monitor) { 87597335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8760dea241dSBarry 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); 87797335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 87897335746SJed Brown } 879cb9d8021SPierre Barbier de Reuille } 880cb9d8021SPierre Barbier de Reuille } else { 8811566a47fSLisandro Dalcin *accept = PETSC_TRUE; 882b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 883cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 884b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 885cb9d8021SPierre Barbier de Reuille } 886cb9d8021SPierre Barbier de Reuille } 887cb9d8021SPierre Barbier de Reuille 8881566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 8891566a47fSLisandro Dalcin PetscReal dt,new_dt; 8901566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 891cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 89297335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 893e6d0a238SBarry Smith adapt->timestepjustincreased += 4; 89497335746SJed Brown if (adapt->monitor) { 89597335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 896e6d0a238SBarry 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); 89797335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 89897335746SJed Brown } 89997335746SJed Brown } 90097335746SJed Brown PetscFunctionReturn(0); 90197335746SJed Brown } 90297335746SJed Brown 90384df9cb4SJed Brown /*@ 90484df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 90584df9cb4SJed Brown 90684df9cb4SJed Brown Collective on MPI_Comm 90784df9cb4SJed Brown 90884df9cb4SJed Brown Input Parameter: 90984df9cb4SJed Brown . comm - The communicator 91084df9cb4SJed Brown 91184df9cb4SJed Brown Output Parameter: 91284df9cb4SJed Brown . adapt - new TSAdapt object 91384df9cb4SJed Brown 91484df9cb4SJed Brown Level: developer 91584df9cb4SJed Brown 91684df9cb4SJed Brown Notes: 91784df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 91884df9cb4SJed Brown 91984df9cb4SJed Brown .keywords: TSAdapt, create 920552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 92184df9cb4SJed Brown @*/ 92284df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 92384df9cb4SJed Brown { 92484df9cb4SJed Brown PetscErrorCode ierr; 92584df9cb4SJed Brown TSAdapt adapt; 92684df9cb4SJed Brown 92784df9cb4SJed Brown PetscFunctionBegin; 9283b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 9293b3bcf4cSLisandro Dalcin *inadapt = NULL; 9303b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 9313b3bcf4cSLisandro Dalcin 93273107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 9331c3436cfSJed Brown 934bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 9351917a363SLisandro Dalcin adapt->safety = 0.9; 9361917a363SLisandro Dalcin adapt->reject_safety = 0.5; 9371917a363SLisandro Dalcin adapt->clip[0] = 0.1; 9381917a363SLisandro Dalcin adapt->clip[1] = 10.; 9391c3436cfSJed Brown adapt->dt_min = 1e-20; 9401917a363SLisandro Dalcin adapt->dt_max = 1e+20; 94197335746SJed Brown adapt->scale_solve_failed = 0.25; 9427619abb3SShri adapt->wnormtype = NORM_2; 9431917a363SLisandro Dalcin 94484df9cb4SJed Brown *inadapt = adapt; 94584df9cb4SJed Brown PetscFunctionReturn(0); 94684df9cb4SJed Brown } 947