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); 128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 13*1917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 1484df9cb4SJed Brown 1584df9cb4SJed Brown /*@C 161c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation 171c84c290SBarry Smith 181c84c290SBarry Smith Not Collective 191c84c290SBarry Smith 201c84c290SBarry Smith Input Parameters: 211c84c290SBarry Smith + name_scheme - name of user-defined adaptivity scheme 221c84c290SBarry Smith - routine_create - routine to create method context 231c84c290SBarry Smith 241c84c290SBarry Smith Notes: 251c84c290SBarry Smith TSAdaptRegister() may be called multiple times to add several user-defined families. 261c84c290SBarry Smith 271c84c290SBarry Smith Sample usage: 281c84c290SBarry Smith .vb 29bdf89e91SBarry Smith TSAdaptRegister("my_scheme",MySchemeCreate); 301c84c290SBarry Smith .ve 311c84c290SBarry Smith 321c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via 331c84c290SBarry Smith $ TSAdaptSetType(ts,"my_scheme") 341c84c290SBarry Smith or at runtime via the option 351c84c290SBarry Smith $ -ts_adapt_type my_scheme 3684df9cb4SJed Brown 3784df9cb4SJed Brown Level: advanced 381c84c290SBarry Smith 391c84c290SBarry Smith .keywords: TSAdapt, register 401c84c290SBarry Smith 411c84c290SBarry Smith .seealso: TSAdaptRegisterAll() 4284df9cb4SJed Brown @*/ 43bdf89e91SBarry Smith PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 4484df9cb4SJed Brown { 4584df9cb4SJed Brown PetscErrorCode ierr; 4684df9cb4SJed Brown 4784df9cb4SJed Brown PetscFunctionBegin; 48a240a19fSJed Brown ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 4984df9cb4SJed Brown PetscFunctionReturn(0); 5084df9cb4SJed Brown } 5184df9cb4SJed Brown 5284df9cb4SJed Brown /*@C 5384df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5484df9cb4SJed Brown 5584df9cb4SJed Brown Not Collective 5684df9cb4SJed Brown 5784df9cb4SJed Brown Level: advanced 5884df9cb4SJed Brown 5984df9cb4SJed Brown .keywords: TSAdapt, register, all 6084df9cb4SJed Brown 6184df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6284df9cb4SJed Brown @*/ 63607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6484df9cb4SJed Brown { 6584df9cb4SJed Brown PetscErrorCode ierr; 6684df9cb4SJed Brown 6784df9cb4SJed Brown PetscFunctionBegin; 680f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 690f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 70bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 711566a47fSLisandro Dalcin ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 72bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 73*1917a363SLisandro Dalcin ierr = TSAdaptRegister(TSADAPTGLEE ,TSAdaptCreate_GLEE);CHKERRQ(ierr); 7484df9cb4SJed Brown PetscFunctionReturn(0); 7584df9cb4SJed Brown } 7684df9cb4SJed Brown 7784df9cb4SJed Brown /*@C 78560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 7984df9cb4SJed Brown called from PetscFinalize(). 8084df9cb4SJed Brown 8184df9cb4SJed Brown Level: developer 8284df9cb4SJed Brown 8384df9cb4SJed Brown .keywords: Petsc, destroy, package 8484df9cb4SJed Brown .seealso: PetscFinalize() 8584df9cb4SJed Brown @*/ 8684df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8784df9cb4SJed Brown { 8837e93019SBarry Smith PetscErrorCode ierr; 8937e93019SBarry Smith 9084df9cb4SJed Brown PetscFunctionBegin; 9137e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9284df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9384df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9484df9cb4SJed Brown PetscFunctionReturn(0); 9584df9cb4SJed Brown } 9684df9cb4SJed Brown 9784df9cb4SJed Brown /*@C 9884df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 9984df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 100*1917a363SLisandro Dalcin TSAdaptCreate() when using static libraries. 10184df9cb4SJed Brown 10284df9cb4SJed Brown Level: developer 10384df9cb4SJed Brown 10484df9cb4SJed Brown .keywords: TSAdapt, initialize, package 10584df9cb4SJed Brown .seealso: PetscInitialize() 10684df9cb4SJed Brown @*/ 107607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 10884df9cb4SJed Brown { 10984df9cb4SJed Brown PetscErrorCode ierr; 11084df9cb4SJed Brown 11184df9cb4SJed Brown PetscFunctionBegin; 11284df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11384df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11484df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 115607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 11684df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 11784df9cb4SJed Brown PetscFunctionReturn(0); 11884df9cb4SJed Brown } 11984df9cb4SJed Brown 1207eef6055SBarry Smith /*@C 1217eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1227eef6055SBarry Smith 1237eef6055SBarry Smith Logicially Collective on TSAdapt 1247eef6055SBarry Smith 1257eef6055SBarry Smith Input Parameter: 126d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt() 1277eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1287eef6055SBarry Smith 1297eef6055SBarry Smith Options Database: 1307eef6055SBarry Smith . -ts_adapt_type basic or none - to setting the adapter type 1317eef6055SBarry Smith 1327eef6055SBarry Smith Level: intermediate 1337eef6055SBarry Smith 1347eef6055SBarry Smith .keywords: TSAdapt, create 1357eef6055SBarry Smith 1367eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1377eef6055SBarry Smith @*/ 13819fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 13984df9cb4SJed Brown { 140ef749922SLisandro Dalcin PetscBool match; 14184df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 14284df9cb4SJed Brown 14384df9cb4SJed Brown PetscFunctionBegin; 1444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 145ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 146ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1471c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 14884df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 149ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1504cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 15184df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15268ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15384df9cb4SJed Brown PetscFunctionReturn(0); 15484df9cb4SJed Brown } 15584df9cb4SJed Brown 156d0288e4fSLisandro Dalcin /*@C 157d0288e4fSLisandro Dalcin TSAdaptGetType - gets the TS adapter method type (as a string). 158d0288e4fSLisandro Dalcin 159d0288e4fSLisandro Dalcin Not Collective 160d0288e4fSLisandro Dalcin 161d0288e4fSLisandro Dalcin Input Parameter: 162d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt() 163d0288e4fSLisandro Dalcin 164d0288e4fSLisandro Dalcin Output Parameter: 165d0288e4fSLisandro Dalcin . type - The name of TS adapter method 166d0288e4fSLisandro Dalcin 167d0288e4fSLisandro Dalcin Level: intermediate 168d0288e4fSLisandro Dalcin 169d0288e4fSLisandro Dalcin .keywords: TSAdapt, get, type 170d0288e4fSLisandro Dalcin .seealso TSAdaptSetType() 171d0288e4fSLisandro Dalcin @*/ 172d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 173d0288e4fSLisandro Dalcin { 174d0288e4fSLisandro Dalcin PetscFunctionBegin; 175d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 176d0288e4fSLisandro Dalcin PetscValidPointer(type,2); 177d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name; 178d0288e4fSLisandro Dalcin PetscFunctionReturn(0); 179d0288e4fSLisandro Dalcin } 180d0288e4fSLisandro Dalcin 18184df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 18284df9cb4SJed Brown { 18384df9cb4SJed Brown PetscErrorCode ierr; 18484df9cb4SJed Brown 18584df9cb4SJed Brown PetscFunctionBegin; 1864782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 18784df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 18884df9cb4SJed Brown PetscFunctionReturn(0); 18984df9cb4SJed Brown } 19084df9cb4SJed Brown 191ad6bc421SBarry Smith /*@C 192ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 193ad6bc421SBarry Smith 194ad6bc421SBarry Smith Collective on PetscViewer 195ad6bc421SBarry Smith 196ad6bc421SBarry Smith Input Parameters: 197ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 198ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 199ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 200ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 201ad6bc421SBarry Smith 202ad6bc421SBarry Smith Level: intermediate 203ad6bc421SBarry Smith 204ad6bc421SBarry Smith Notes: 205ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 206ad6bc421SBarry Smith 207ad6bc421SBarry Smith Notes for advanced users: 208ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 209ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 210ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 211ad6bc421SBarry Smith format is 212ad6bc421SBarry Smith .vb 213ad6bc421SBarry Smith has not yet been determined 214ad6bc421SBarry Smith .ve 215ad6bc421SBarry Smith 216ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 217ad6bc421SBarry Smith @*/ 2184782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 219ad6bc421SBarry Smith { 220ad6bc421SBarry Smith PetscErrorCode ierr; 221ad6bc421SBarry Smith PetscBool isbinary; 222ad6bc421SBarry Smith char type[256]; 223ad6bc421SBarry Smith 224ad6bc421SBarry Smith PetscFunctionBegin; 2254782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 226ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 227ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 228ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 229ad6bc421SBarry Smith 230060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2314782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2324782b174SLisandro Dalcin if (adapt->ops->load) { 2334782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 234ad6bc421SBarry Smith } 235ad6bc421SBarry Smith PetscFunctionReturn(0); 236ad6bc421SBarry Smith } 237ad6bc421SBarry Smith 23884df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 23984df9cb4SJed Brown { 24084df9cb4SJed Brown PetscErrorCode ierr; 241*1917a363SLisandro Dalcin PetscBool iascii,isbinary,isnone; 24284df9cb4SJed Brown 24384df9cb4SJed Brown PetscFunctionBegin; 2444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2454782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2464782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2474782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 248251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 249ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 25084df9cb4SJed Brown if (iascii) { 251dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 252*1917a363SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr); 253*1917a363SLisandro Dalcin if (!isnone) { 254bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 255*1917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 256*1917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 257*1917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 258*1917a363SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 259bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 260bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 261*1917a363SLisandro Dalcin } 26284df9cb4SJed Brown if (adapt->ops->view) { 26384df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 26484df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 26584df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 26684df9cb4SJed Brown } 267ad6bc421SBarry Smith } else if (isbinary) { 268ad6bc421SBarry Smith char type[256]; 269ad6bc421SBarry Smith 270ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 271ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 272ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 273bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 274f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 275f2c2a1b9SBarry Smith } 27684df9cb4SJed Brown PetscFunctionReturn(0); 27784df9cb4SJed Brown } 27884df9cb4SJed Brown 279099af0a3SLisandro Dalcin /*@ 280099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 281099af0a3SLisandro Dalcin 282099af0a3SLisandro Dalcin Collective on TS 283099af0a3SLisandro Dalcin 284099af0a3SLisandro Dalcin Input Parameter: 285099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 286099af0a3SLisandro Dalcin 287099af0a3SLisandro Dalcin Level: developer 288099af0a3SLisandro Dalcin 289099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 290099af0a3SLisandro Dalcin @*/ 291099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 292099af0a3SLisandro Dalcin { 293099af0a3SLisandro Dalcin PetscErrorCode ierr; 294099af0a3SLisandro Dalcin 295099af0a3SLisandro Dalcin PetscFunctionBegin; 296099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 297099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 298099af0a3SLisandro Dalcin PetscFunctionReturn(0); 299099af0a3SLisandro Dalcin } 300099af0a3SLisandro Dalcin 30184df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 30284df9cb4SJed Brown { 30384df9cb4SJed Brown PetscErrorCode ierr; 30484df9cb4SJed Brown 30584df9cb4SJed Brown PetscFunctionBegin; 30684df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 30784df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 3084782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 309099af0a3SLisandro Dalcin 310099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 311099af0a3SLisandro Dalcin 31284df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 3131c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 31484df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 31584df9cb4SJed Brown PetscFunctionReturn(0); 31684df9cb4SJed Brown } 31784df9cb4SJed Brown 3181c3436cfSJed Brown /*@ 3191c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 3201c3436cfSJed Brown 3211c3436cfSJed Brown Collective on TSAdapt 3221c3436cfSJed Brown 3231c3436cfSJed Brown Input Arguments: 3241c3436cfSJed Brown + adapt - adaptive controller context 3251c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 3261c3436cfSJed Brown 327bf997491SLisandro Dalcin Options Database Keys: 328bf997491SLisandro Dalcin . -ts_adapt_monitor 329bf997491SLisandro Dalcin 3301c3436cfSJed Brown Level: intermediate 3311c3436cfSJed Brown 3321c3436cfSJed Brown .seealso: TSAdaptChoose() 3331c3436cfSJed Brown @*/ 3341c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3351c3436cfSJed Brown { 3361c3436cfSJed Brown PetscErrorCode ierr; 3371c3436cfSJed Brown 3381c3436cfSJed Brown PetscFunctionBegin; 3394782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3404782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3411c3436cfSJed Brown if (flg) { 342ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3431c3436cfSJed Brown } else { 3441c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3451c3436cfSJed Brown } 3461c3436cfSJed Brown PetscFunctionReturn(0); 3471c3436cfSJed Brown } 3481c3436cfSJed Brown 3490873808bSJed Brown /*@C 350bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3510873808bSJed Brown 3520873808bSJed Brown Logically collective on TSAdapt 3530873808bSJed Brown 3540873808bSJed Brown Input Arguments: 3550873808bSJed Brown + adapt - adaptive controller context 3560873808bSJed Brown - func - stage check function 3570873808bSJed Brown 3580873808bSJed Brown Arguments of func: 3590873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3600873808bSJed Brown 3610873808bSJed Brown + adapt - adaptive controller context 3620873808bSJed Brown . ts - time stepping context 3630873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3640873808bSJed Brown 3650873808bSJed Brown Level: advanced 3660873808bSJed Brown 3670873808bSJed Brown .seealso: TSAdaptChoose() 3680873808bSJed Brown @*/ 369b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3700873808bSJed Brown { 3710873808bSJed Brown 3720873808bSJed Brown PetscFunctionBegin; 3730873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 37468ae941aSLisandro Dalcin adapt->checkstage = func; 3750873808bSJed Brown PetscFunctionReturn(0); 3760873808bSJed Brown } 3770873808bSJed Brown 3781c3436cfSJed Brown /*@ 379bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 380bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 381bf997491SLisandro Dalcin 382*1917a363SLisandro Dalcin Logically collective on TSAdapt 383bf997491SLisandro Dalcin 384bf997491SLisandro Dalcin Input Arguments: 385bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 386bf997491SLisandro Dalcin - flag - whether to always accept steps 387bf997491SLisandro Dalcin 388bf997491SLisandro Dalcin Options Database Keys: 389bf997491SLisandro Dalcin . -ts_adapt_always_accept 390bf997491SLisandro Dalcin 391bf997491SLisandro Dalcin Level: intermediate 392bf997491SLisandro Dalcin 393bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 394bf997491SLisandro Dalcin @*/ 395bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 396bf997491SLisandro Dalcin { 397bf997491SLisandro Dalcin PetscFunctionBegin; 398bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 399bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 400bf997491SLisandro Dalcin adapt->always_accept = flag; 401bf997491SLisandro Dalcin PetscFunctionReturn(0); 402bf997491SLisandro Dalcin } 403bf997491SLisandro Dalcin 404bf997491SLisandro Dalcin /*@ 405*1917a363SLisandro Dalcin TSAdaptSetSafety - Set safety factors 4061c3436cfSJed Brown 407*1917a363SLisandro Dalcin Logically collective on TSAdapt 408*1917a363SLisandro Dalcin 409*1917a363SLisandro Dalcin Input Arguments: 410*1917a363SLisandro Dalcin + adapt - adaptive controller context 411*1917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 412*1917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 413*1917a363SLisandro Dalcin 414*1917a363SLisandro Dalcin Options Database Keys: 415*1917a363SLisandro Dalcin + -ts_adapt_safety 416*1917a363SLisandro Dalcin - -ts_adapt_reject_safety 417*1917a363SLisandro Dalcin 418*1917a363SLisandro Dalcin Level: intermediate 419*1917a363SLisandro Dalcin 420*1917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 421*1917a363SLisandro Dalcin @*/ 422*1917a363SLisandro Dalcin PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 423*1917a363SLisandro Dalcin { 424*1917a363SLisandro Dalcin PetscFunctionBegin; 425*1917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 426*1917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,safety,2); 427*1917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 428*1917a363SLisandro Dalcin if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 429*1917a363SLisandro 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); 430*1917a363SLisandro 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); 431*1917a363SLisandro 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); 432*1917a363SLisandro Dalcin if (safety != PETSC_DEFAULT) adapt->safety = safety; 433*1917a363SLisandro Dalcin if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 434*1917a363SLisandro Dalcin PetscFunctionReturn(0); 435*1917a363SLisandro Dalcin } 436*1917a363SLisandro Dalcin 437*1917a363SLisandro Dalcin /*@ 438*1917a363SLisandro Dalcin TSAdaptGetSafety - Get safety factors 439*1917a363SLisandro Dalcin 440*1917a363SLisandro Dalcin Not Collective 441*1917a363SLisandro Dalcin 442*1917a363SLisandro Dalcin Input Arguments: 443*1917a363SLisandro Dalcin . adapt - adaptive controller context 444*1917a363SLisandro Dalcin 445*1917a363SLisandro Dalcin Ouput Arguments: 446*1917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal 447*1917a363SLisandro Dalcin + reject_safety - extra safety factor to apply if the last step was rejected 448*1917a363SLisandro Dalcin 449*1917a363SLisandro Dalcin Level: intermediate 450*1917a363SLisandro Dalcin 451*1917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 452*1917a363SLisandro Dalcin @*/ 453*1917a363SLisandro Dalcin PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 454*1917a363SLisandro Dalcin { 455*1917a363SLisandro Dalcin PetscFunctionBegin; 456*1917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 457*1917a363SLisandro Dalcin if (safety) PetscValidRealPointer(safety,2); 458*1917a363SLisandro Dalcin if (reject_safety) PetscValidRealPointer(reject_safety,3); 459*1917a363SLisandro Dalcin if (safety) *safety = adapt->safety; 460*1917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety; 461*1917a363SLisandro Dalcin PetscFunctionReturn(0); 462*1917a363SLisandro Dalcin } 463*1917a363SLisandro Dalcin 464*1917a363SLisandro Dalcin /*@ 465*1917a363SLisandro Dalcin TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 466*1917a363SLisandro Dalcin 467*1917a363SLisandro Dalcin Logically collective on TSAdapt 468*1917a363SLisandro Dalcin 469*1917a363SLisandro Dalcin Input Arguments: 470*1917a363SLisandro Dalcin + adapt - adaptive controller context 471*1917a363SLisandro Dalcin . low - admissible decrease factor 472*1917a363SLisandro Dalcin + high - admissible increase factor 473*1917a363SLisandro Dalcin 474*1917a363SLisandro Dalcin Options Database Keys: 475*1917a363SLisandro Dalcin . -ts_adapt_clip 476*1917a363SLisandro Dalcin 477*1917a363SLisandro Dalcin Level: intermediate 478*1917a363SLisandro Dalcin 479*1917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptGetClip() 480*1917a363SLisandro Dalcin @*/ 481*1917a363SLisandro Dalcin PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 482*1917a363SLisandro Dalcin { 483*1917a363SLisandro Dalcin PetscFunctionBegin; 484*1917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 485*1917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,low,2); 486*1917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,high,3); 487*1917a363SLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 488*1917a363SLisandro 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); 489*1917a363SLisandro 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); 490*1917a363SLisandro Dalcin if (low != PETSC_DEFAULT) adapt->clip[0] = low; 491*1917a363SLisandro Dalcin if (high != PETSC_DEFAULT) adapt->clip[1] = high; 492*1917a363SLisandro Dalcin PetscFunctionReturn(0); 493*1917a363SLisandro Dalcin } 494*1917a363SLisandro Dalcin 495*1917a363SLisandro Dalcin /*@ 496*1917a363SLisandro Dalcin TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 497*1917a363SLisandro Dalcin 498*1917a363SLisandro Dalcin Not Collective 499*1917a363SLisandro Dalcin 500*1917a363SLisandro Dalcin Input Arguments: 501*1917a363SLisandro Dalcin . adapt - adaptive controller context 502*1917a363SLisandro Dalcin 503*1917a363SLisandro Dalcin Ouput Arguments: 504*1917a363SLisandro Dalcin + low - optional, admissible decrease factor 505*1917a363SLisandro Dalcin - high - optional, admissible increase factor 506*1917a363SLisandro Dalcin 507*1917a363SLisandro Dalcin Level: intermediate 508*1917a363SLisandro Dalcin 509*1917a363SLisandro Dalcin .seealso: TSAdaptChoose(), TSAdaptSetClip() 510*1917a363SLisandro Dalcin @*/ 511*1917a363SLisandro Dalcin PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 512*1917a363SLisandro Dalcin { 513*1917a363SLisandro Dalcin PetscFunctionBegin; 514*1917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 515*1917a363SLisandro Dalcin if (low) PetscValidRealPointer(low,2); 516*1917a363SLisandro Dalcin if (high) PetscValidRealPointer(high,3); 517*1917a363SLisandro Dalcin if (low) *low = adapt->clip[0]; 518*1917a363SLisandro Dalcin if (high) *high = adapt->clip[1]; 519*1917a363SLisandro Dalcin PetscFunctionReturn(0); 520*1917a363SLisandro Dalcin } 521*1917a363SLisandro Dalcin 522*1917a363SLisandro Dalcin /*@ 523*1917a363SLisandro Dalcin TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 524*1917a363SLisandro Dalcin 525*1917a363SLisandro Dalcin Logically collective on TSAdapt 5261c3436cfSJed Brown 5271c3436cfSJed Brown Input Arguments: 528552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 5291c3436cfSJed Brown . hmin - minimum time step 5301c3436cfSJed Brown - hmax - maximum time step 5311c3436cfSJed Brown 5321c3436cfSJed Brown Options Database Keys: 5331c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 5341c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 5351c3436cfSJed Brown 5361c3436cfSJed Brown Level: intermediate 5371c3436cfSJed Brown 538*1917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 5391c3436cfSJed Brown @*/ 5401c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 5411c3436cfSJed Brown { 5421c3436cfSJed Brown 5431c3436cfSJed Brown PetscFunctionBegin; 5444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 545b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 546b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 547b1f5bccaSLisandro 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); 548b1f5bccaSLisandro 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); 549b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 550b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 551b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 552b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 553b1f5bccaSLisandro 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); 554*1917a363SLisandro Dalcin PetscFunctionReturn(0); 555*1917a363SLisandro Dalcin } 556*1917a363SLisandro Dalcin 557*1917a363SLisandro Dalcin /*@ 558*1917a363SLisandro Dalcin TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 559*1917a363SLisandro Dalcin 560*1917a363SLisandro Dalcin Not Collective 561*1917a363SLisandro Dalcin 562*1917a363SLisandro Dalcin Input Arguments: 563*1917a363SLisandro Dalcin . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 564*1917a363SLisandro Dalcin 565*1917a363SLisandro Dalcin Output Arguments: 566*1917a363SLisandro Dalcin + hmin - minimum time step 567*1917a363SLisandro Dalcin - hmax - maximum time step 568*1917a363SLisandro Dalcin 569*1917a363SLisandro Dalcin Level: intermediate 570*1917a363SLisandro Dalcin 571*1917a363SLisandro Dalcin .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 572*1917a363SLisandro Dalcin @*/ 573*1917a363SLisandro Dalcin PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 574*1917a363SLisandro Dalcin { 575*1917a363SLisandro Dalcin 576*1917a363SLisandro Dalcin PetscFunctionBegin; 577*1917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 578*1917a363SLisandro Dalcin if (hmin) PetscValidRealPointer(hmin,2); 579*1917a363SLisandro Dalcin if (hmax) PetscValidRealPointer(hmax,3); 580*1917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min; 581*1917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max; 5821c3436cfSJed Brown PetscFunctionReturn(0); 5831c3436cfSJed Brown } 5841c3436cfSJed Brown 585e55864a3SBarry Smith /* 58684df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 58784df9cb4SJed Brown 58884df9cb4SJed Brown Collective on TSAdapt 58984df9cb4SJed Brown 59084df9cb4SJed Brown Input Parameter: 59184df9cb4SJed Brown . adapt - the TSAdapt context 59284df9cb4SJed Brown 59384df9cb4SJed Brown Options Database Keys: 594*1917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity 595bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 596*1917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 597*1917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 598*1917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 59923de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 60023de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 60123de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 60223de1d84SBarry Smith - -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 60384df9cb4SJed Brown 60484df9cb4SJed Brown Level: advanced 60584df9cb4SJed Brown 60684df9cb4SJed Brown Notes: 60784df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 60884df9cb4SJed Brown 60923de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 61084df9cb4SJed Brown 611*1917a363SLisandro Dalcin .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 612*1917a363SLisandro Dalcin TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 613e55864a3SBarry Smith */ 6144416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 61584df9cb4SJed Brown { 61684df9cb4SJed Brown PetscErrorCode ierr; 61784df9cb4SJed Brown char type[256] = TSADAPTBASIC; 618*1917a363SLisandro Dalcin PetscReal safety,reject_safety,clip[2],hmin,hmax; 6191c3436cfSJed Brown PetscBool set,flg; 620*1917a363SLisandro Dalcin PetscInt two; 62184df9cb4SJed Brown 62284df9cb4SJed Brown PetscFunctionBegin; 6234782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 62484df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 6251566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 626e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 62723de1d84SBarry 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); 62884df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 62984df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 63084df9cb4SJed Brown } 631*1917a363SLisandro Dalcin 632bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 633bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 634*1917a363SLisandro Dalcin 635*1917a363SLisandro Dalcin safety = adapt->safety; reject_safety = adapt->reject_safety; 636*1917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 637*1917a363SLisandro 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); 638*1917a363SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 639*1917a363SLisandro Dalcin 640*1917a363SLisandro Dalcin two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 641*1917a363SLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 642*1917a363SLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 643*1917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 644*1917a363SLisandro Dalcin 645*1917a363SLisandro Dalcin hmin = adapt->dt_min; hmax = adapt->dt_max; 646*1917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 647*1917a363SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 648bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 649*1917a363SLisandro Dalcin 6500298fd71SBarry 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); 651*1917a363SLisandro Dalcin 6527619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 6537619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 654*1917a363SLisandro Dalcin 655*1917a363SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 656*1917a363SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 657*1917a363SLisandro Dalcin 658e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 65984df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 66084df9cb4SJed Brown PetscFunctionReturn(0); 66184df9cb4SJed Brown } 66284df9cb4SJed Brown 66384df9cb4SJed Brown /*@ 66484df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 66584df9cb4SJed Brown 666*1917a363SLisandro Dalcin Logically collective on TSAdapt 66784df9cb4SJed Brown 66884df9cb4SJed Brown Input Argument: 66984df9cb4SJed Brown . adapt - adaptive controller 67084df9cb4SJed Brown 67184df9cb4SJed Brown Level: developer 67284df9cb4SJed Brown 67384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 67484df9cb4SJed Brown @*/ 67584df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 67684df9cb4SJed Brown { 67784df9cb4SJed Brown PetscErrorCode ierr; 67884df9cb4SJed Brown 67984df9cb4SJed Brown PetscFunctionBegin; 6804782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 68184df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 68284df9cb4SJed Brown PetscFunctionReturn(0); 68384df9cb4SJed Brown } 68484df9cb4SJed Brown 68584df9cb4SJed Brown /*@C 68684df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 68784df9cb4SJed Brown 688*1917a363SLisandro Dalcin Logically collective on TSAdapt 68984df9cb4SJed Brown 69084df9cb4SJed Brown Input Arguments: 691552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 69284df9cb4SJed Brown . name - name of the candidate scheme to add 69384df9cb4SJed Brown . order - order of the candidate scheme 69484df9cb4SJed Brown . stageorder - stage order of the candidate scheme 6958d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 69684df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 69784df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 69884df9cb4SJed Brown 69984df9cb4SJed Brown Note: 70084df9cb4SJed Brown This routine is not available in Fortran. 70184df9cb4SJed Brown 70284df9cb4SJed Brown Level: developer 70384df9cb4SJed Brown 70484df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 70584df9cb4SJed Brown @*/ 7068d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 70784df9cb4SJed Brown { 70884df9cb4SJed Brown PetscInt c; 70984df9cb4SJed Brown 71084df9cb4SJed Brown PetscFunctionBegin; 71184df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 712ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 71384df9cb4SJed Brown if (inuse) { 714ce94432eSBarry 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()"); 71584df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 71684df9cb4SJed Brown } 7171c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 7181c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 719bbd56ea5SKarl Rupp 72084df9cb4SJed Brown adapt->candidates.name[c] = name; 72184df9cb4SJed Brown adapt->candidates.order[c] = order; 72284df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 7238d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 72484df9cb4SJed Brown adapt->candidates.cost[c] = cost; 72584df9cb4SJed Brown adapt->candidates.n++; 72684df9cb4SJed Brown PetscFunctionReturn(0); 72784df9cb4SJed Brown } 72884df9cb4SJed Brown 7298d59e960SJed Brown /*@C 7308d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 7318d59e960SJed Brown 7328d59e960SJed Brown Not Collective 7338d59e960SJed Brown 7348d59e960SJed Brown Input Arguments: 7358d59e960SJed Brown . adapt - time step adaptivity context 7368d59e960SJed Brown 7378d59e960SJed Brown Output Arguments: 7388d59e960SJed Brown + n - number of candidate schemes, always at least 1 7398d59e960SJed Brown . order - the order of each candidate scheme 7408d59e960SJed Brown . stageorder - the stage order of each candidate scheme 7418d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 7428d59e960SJed Brown - cost - the relative cost of each scheme 7438d59e960SJed Brown 7448d59e960SJed Brown Level: developer 7458d59e960SJed Brown 7468d59e960SJed Brown Note: 7478d59e960SJed Brown The current scheme is always returned in the first slot 7488d59e960SJed Brown 7498d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 7508d59e960SJed Brown @*/ 7518d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 7528d59e960SJed Brown { 7538d59e960SJed Brown PetscFunctionBegin; 7548d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 7558d59e960SJed Brown if (n) *n = adapt->candidates.n; 7568d59e960SJed Brown if (order) *order = adapt->candidates.order; 7578d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 7588d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 7598d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 7608d59e960SJed Brown PetscFunctionReturn(0); 7618d59e960SJed Brown } 7628d59e960SJed Brown 76384df9cb4SJed Brown /*@C 76484df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 76584df9cb4SJed Brown 766*1917a363SLisandro Dalcin Collective on TSAdapt 76784df9cb4SJed Brown 76884df9cb4SJed Brown Input Arguments: 76984df9cb4SJed Brown + adapt - adaptive contoller 77084df9cb4SJed Brown - h - current step size 77184df9cb4SJed Brown 77284df9cb4SJed Brown Output Arguments: 7731566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 77484df9cb4SJed Brown . next_h - step size to use for the next step 77584df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 77684df9cb4SJed Brown 7771c3436cfSJed Brown Note: 7781c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 7791c3436cfSJed Brown being retried after an initial rejection. 7801c3436cfSJed Brown 78184df9cb4SJed Brown Level: developer 78284df9cb4SJed Brown 78384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 78484df9cb4SJed Brown @*/ 78584df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 78684df9cb4SJed Brown { 78784df9cb4SJed Brown PetscErrorCode ierr; 7881566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 7891566a47fSLisandro Dalcin PetscInt scheme = 0; 7900b99f514SJed Brown PetscReal wlte = -1.0; 7915e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 7925e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 79384df9cb4SJed Brown 79484df9cb4SJed Brown PetscFunctionBegin; 79584df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 79684df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 7971566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 79884df9cb4SJed Brown PetscValidPointer(next_h,5); 79984df9cb4SJed Brown PetscValidIntPointer(accept,6); 8001566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 8011566a47fSLisandro Dalcin 8021566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 8031566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 8041566a47fSLisandro Dalcin *next_h = h; 8051566a47fSLisandro Dalcin *accept = PETSC_TRUE; 8061566a47fSLisandro Dalcin PetscFunctionReturn(0); 8071566a47fSLisandro Dalcin } 8081566a47fSLisandro Dalcin 8095e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 8101566a47fSLisandro 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); 8111566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 8121566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 8131566a47fSLisandro Dalcin 81449354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 81549354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 8161566a47fSLisandro Dalcin if (ts->ptime + ts->time_step + *next_h >= ts->max_time) { 8171566a47fSLisandro Dalcin PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step); 8188709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 81949354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 82049354f04SShri Abhyankar } 82149354f04SShri Abhyankar } 8221c3436cfSJed Brown 8231c3436cfSJed Brown if (adapt->monitor) { 8241566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 8251c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8260b99f514SJed Brown if (wlte < 0) { 8270dea241dSBarry 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); 8280b99f514SJed Brown } else { 8290dea241dSBarry 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); 8300b99f514SJed Brown } 8311c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8321c3436cfSJed Brown } 83384df9cb4SJed Brown PetscFunctionReturn(0); 83484df9cb4SJed Brown } 83584df9cb4SJed Brown 83697335746SJed Brown /*@ 83797335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 83897335746SJed Brown 839*1917a363SLisandro Dalcin Collective on TSAdapt 84097335746SJed Brown 84197335746SJed Brown Input Arguments: 84297335746SJed Brown + adapt - adaptive controller context 843b295832fSPierre Barbier de Reuille . ts - time stepper 844b295832fSPierre Barbier de Reuille . t - Current simulation time 845b295832fSPierre Barbier de Reuille - Y - Current solution vector 84697335746SJed Brown 84797335746SJed Brown Output Arguments: 84897335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 84997335746SJed Brown 85097335746SJed Brown Level: developer 85197335746SJed Brown 85297335746SJed Brown .seealso: 85397335746SJed Brown @*/ 854b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 85597335746SJed Brown { 85697335746SJed Brown PetscErrorCode ierr; 8571566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 85897335746SJed Brown 85997335746SJed Brown PetscFunctionBegin; 8604782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 8614782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 8624782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 8631566a47fSLisandro Dalcin 8641566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 86597335746SJed Brown if (snesreason < 0) { 86697335746SJed Brown *accept = PETSC_FALSE; 8676de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 86897335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 86948c19aefSLisandro 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); 87097335746SJed Brown if (adapt->monitor) { 87197335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8720dea241dSBarry 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); 87397335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 87497335746SJed Brown } 875cb9d8021SPierre Barbier de Reuille } 876cb9d8021SPierre Barbier de Reuille } else { 8771566a47fSLisandro Dalcin *accept = PETSC_TRUE; 878b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 879cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 880b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 881cb9d8021SPierre Barbier de Reuille } 882cb9d8021SPierre Barbier de Reuille } 883cb9d8021SPierre Barbier de Reuille 8841566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 8851566a47fSLisandro Dalcin PetscReal dt,new_dt; 8861566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 887cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 88897335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 88997335746SJed Brown if (adapt->monitor) { 89097335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 8910dea241dSBarry Smith ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr); 89297335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 89397335746SJed Brown } 89497335746SJed Brown } 89597335746SJed Brown PetscFunctionReturn(0); 89697335746SJed Brown } 89797335746SJed Brown 89884df9cb4SJed Brown /*@ 89984df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 90084df9cb4SJed Brown 90184df9cb4SJed Brown Collective on MPI_Comm 90284df9cb4SJed Brown 90384df9cb4SJed Brown Input Parameter: 90484df9cb4SJed Brown . comm - The communicator 90584df9cb4SJed Brown 90684df9cb4SJed Brown Output Parameter: 90784df9cb4SJed Brown . adapt - new TSAdapt object 90884df9cb4SJed Brown 90984df9cb4SJed Brown Level: developer 91084df9cb4SJed Brown 91184df9cb4SJed Brown Notes: 91284df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 91384df9cb4SJed Brown 91484df9cb4SJed Brown .keywords: TSAdapt, create 915552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 91684df9cb4SJed Brown @*/ 91784df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 91884df9cb4SJed Brown { 91984df9cb4SJed Brown PetscErrorCode ierr; 92084df9cb4SJed Brown TSAdapt adapt; 92184df9cb4SJed Brown 92284df9cb4SJed Brown PetscFunctionBegin; 9233b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 9243b3bcf4cSLisandro Dalcin *inadapt = NULL; 9253b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 9263b3bcf4cSLisandro Dalcin 92773107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 9281c3436cfSJed Brown 929bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 930*1917a363SLisandro Dalcin adapt->safety = 0.9; 931*1917a363SLisandro Dalcin adapt->reject_safety = 0.5; 932*1917a363SLisandro Dalcin adapt->clip[0] = 0.1; 933*1917a363SLisandro Dalcin adapt->clip[1] = 10.; 9341c3436cfSJed Brown adapt->dt_min = 1e-20; 935*1917a363SLisandro Dalcin adapt->dt_max = 1e+20; 93697335746SJed Brown adapt->scale_solve_failed = 0.25; 9377619abb3SShri adapt->wnormtype = NORM_2; 938*1917a363SLisandro Dalcin 9397eef6055SBarry Smith ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 9401c3436cfSJed Brown 94184df9cb4SJed Brown *inadapt = adapt; 94284df9cb4SJed Brown PetscFunctionReturn(0); 94384df9cb4SJed Brown } 944