184df9cb4SJed Brown 2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 384df9cb4SJed Brown 4140e18c1SBarry Smith static PetscFunctionList TSAdaptList; 584df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 684df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 784df9cb4SJed Brown static PetscClassId TSADAPT_CLASSID; 884df9cb4SJed Brown 98cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 1284df9cb4SJed Brown 1384df9cb4SJed Brown #undef __FUNCT__ 1484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister" 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 #undef __FUNCT__ 5384df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll" 5484df9cb4SJed Brown /*@C 5584df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5684df9cb4SJed Brown 5784df9cb4SJed Brown Not Collective 5884df9cb4SJed Brown 5984df9cb4SJed Brown Level: advanced 6084df9cb4SJed Brown 6184df9cb4SJed Brown .keywords: TSAdapt, register, all 6284df9cb4SJed Brown 6384df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6484df9cb4SJed Brown @*/ 65607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6684df9cb4SJed Brown { 6784df9cb4SJed Brown PetscErrorCode ierr; 6884df9cb4SJed Brown 6984df9cb4SJed Brown PetscFunctionBegin; 700f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 710f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 72bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 73bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 74bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 7584df9cb4SJed Brown PetscFunctionReturn(0); 7684df9cb4SJed Brown } 7784df9cb4SJed Brown 7884df9cb4SJed Brown #undef __FUNCT__ 7984df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage" 8084df9cb4SJed Brown /*@C 8184df9cb4SJed Brown TSFinalizePackage - This function destroys everything in the TS package. It is 8284df9cb4SJed Brown called from PetscFinalize(). 8384df9cb4SJed Brown 8484df9cb4SJed Brown Level: developer 8584df9cb4SJed Brown 8684df9cb4SJed Brown .keywords: Petsc, destroy, package 8784df9cb4SJed Brown .seealso: PetscFinalize() 8884df9cb4SJed Brown @*/ 8984df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 9084df9cb4SJed Brown { 9137e93019SBarry Smith PetscErrorCode ierr; 9237e93019SBarry Smith 9384df9cb4SJed Brown PetscFunctionBegin; 9437e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9584df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9684df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9784df9cb4SJed Brown PetscFunctionReturn(0); 9884df9cb4SJed Brown } 9984df9cb4SJed Brown 10084df9cb4SJed Brown #undef __FUNCT__ 10184df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage" 10284df9cb4SJed Brown /*@C 10384df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 10484df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 10584df9cb4SJed Brown TSCreate_GL() when using static libraries. 10684df9cb4SJed Brown 10784df9cb4SJed Brown Level: developer 10884df9cb4SJed Brown 10984df9cb4SJed Brown .keywords: TSAdapt, initialize, package 11084df9cb4SJed Brown .seealso: PetscInitialize() 11184df9cb4SJed Brown @*/ 112607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11384df9cb4SJed Brown { 11484df9cb4SJed Brown PetscErrorCode ierr; 11584df9cb4SJed Brown 11684df9cb4SJed Brown PetscFunctionBegin; 11784df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11884df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11984df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 120607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 12184df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 12284df9cb4SJed Brown PetscFunctionReturn(0); 12384df9cb4SJed Brown } 12484df9cb4SJed Brown 12584df9cb4SJed Brown #undef __FUNCT__ 12684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 12719fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 12884df9cb4SJed Brown { 129ef749922SLisandro Dalcin PetscBool match; 1303403ddbcSLisandro Dalcin PetscErrorCode (*checkstage)(TSAdapt,TS,PetscBool*); 13184df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 13284df9cb4SJed Brown 13384df9cb4SJed Brown PetscFunctionBegin; 1344782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 135ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 136ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1371c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 13884df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 139ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1404cefc2ffSBarry Smith /* Reinitialize function pointers in TSAdaptOps structure */ 1413403ddbcSLisandro Dalcin checkstage = adapt->ops->checkstage; 1424cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 1433403ddbcSLisandro Dalcin adapt->ops->checkstage = checkstage; 14484df9cb4SJed Brown ierr = (*r)(adapt);CHKERRQ(ierr); 14584df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 14684df9cb4SJed Brown PetscFunctionReturn(0); 14784df9cb4SJed Brown } 14884df9cb4SJed Brown 14984df9cb4SJed Brown #undef __FUNCT__ 15084df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 15184df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 15284df9cb4SJed Brown { 15384df9cb4SJed Brown PetscErrorCode ierr; 15484df9cb4SJed Brown 15584df9cb4SJed Brown PetscFunctionBegin; 1564782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 15784df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 15884df9cb4SJed Brown PetscFunctionReturn(0); 15984df9cb4SJed Brown } 16084df9cb4SJed Brown 16184df9cb4SJed Brown #undef __FUNCT__ 162ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad" 163ad6bc421SBarry Smith /*@C 164ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 165ad6bc421SBarry Smith 166ad6bc421SBarry Smith Collective on PetscViewer 167ad6bc421SBarry Smith 168ad6bc421SBarry Smith Input Parameters: 169ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 170ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 171ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 172ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 173ad6bc421SBarry Smith 174ad6bc421SBarry Smith Level: intermediate 175ad6bc421SBarry Smith 176ad6bc421SBarry Smith Notes: 177ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 178ad6bc421SBarry Smith 179ad6bc421SBarry Smith Notes for advanced users: 180ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 181ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 182ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 183ad6bc421SBarry Smith format is 184ad6bc421SBarry Smith .vb 185ad6bc421SBarry Smith has not yet been determined 186ad6bc421SBarry Smith .ve 187ad6bc421SBarry Smith 188ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 189ad6bc421SBarry Smith @*/ 1904782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 191ad6bc421SBarry Smith { 192ad6bc421SBarry Smith PetscErrorCode ierr; 193ad6bc421SBarry Smith PetscBool isbinary; 194ad6bc421SBarry Smith char type[256]; 195ad6bc421SBarry Smith 196ad6bc421SBarry Smith PetscFunctionBegin; 1974782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 198ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 199ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 200ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 201ad6bc421SBarry Smith 202ad6bc421SBarry Smith ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 2034782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt, type);CHKERRQ(ierr); 2044782b174SLisandro Dalcin if (adapt->ops->load) { 2054782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 206ad6bc421SBarry Smith } 207ad6bc421SBarry Smith PetscFunctionReturn(0); 208ad6bc421SBarry Smith } 209ad6bc421SBarry Smith 210ad6bc421SBarry Smith #undef __FUNCT__ 21184df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 21284df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 21384df9cb4SJed Brown { 21484df9cb4SJed Brown PetscErrorCode ierr; 215ad6bc421SBarry Smith PetscBool iascii,isbinary; 21684df9cb4SJed Brown 21784df9cb4SJed Brown PetscFunctionBegin; 2184782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2194782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2204782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2214782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 222251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 223ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 22484df9cb4SJed Brown if (iascii) { 225dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 22684df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 22784df9cb4SJed Brown if (adapt->ops->view) { 22884df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 22984df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 23084df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 23184df9cb4SJed Brown } 232ad6bc421SBarry Smith } else if (isbinary) { 233ad6bc421SBarry Smith char type[256]; 234ad6bc421SBarry Smith 235ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 236ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 237ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 238bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 239f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 240f2c2a1b9SBarry Smith } 24184df9cb4SJed Brown PetscFunctionReturn(0); 24284df9cb4SJed Brown } 24384df9cb4SJed Brown 24484df9cb4SJed Brown #undef __FUNCT__ 245*099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset" 246*099af0a3SLisandro Dalcin /*@ 247*099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 248*099af0a3SLisandro Dalcin 249*099af0a3SLisandro Dalcin Collective on TS 250*099af0a3SLisandro Dalcin 251*099af0a3SLisandro Dalcin Input Parameter: 252*099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 253*099af0a3SLisandro Dalcin 254*099af0a3SLisandro Dalcin Level: developer 255*099af0a3SLisandro Dalcin 256*099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 257*099af0a3SLisandro Dalcin @*/ 258*099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 259*099af0a3SLisandro Dalcin { 260*099af0a3SLisandro Dalcin PetscErrorCode ierr; 261*099af0a3SLisandro Dalcin 262*099af0a3SLisandro Dalcin PetscFunctionBegin; 263*099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 264*099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 265*099af0a3SLisandro Dalcin PetscFunctionReturn(0); 266*099af0a3SLisandro Dalcin } 267*099af0a3SLisandro Dalcin 268*099af0a3SLisandro Dalcin #undef __FUNCT__ 26984df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 27084df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 27184df9cb4SJed Brown { 27284df9cb4SJed Brown PetscErrorCode ierr; 27384df9cb4SJed Brown 27484df9cb4SJed Brown PetscFunctionBegin; 27584df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 27684df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 2774782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 278*099af0a3SLisandro Dalcin 279*099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 280*099af0a3SLisandro Dalcin 28184df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 2821c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 28384df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 28484df9cb4SJed Brown PetscFunctionReturn(0); 28584df9cb4SJed Brown } 28684df9cb4SJed Brown 28784df9cb4SJed Brown #undef __FUNCT__ 2881c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 2891c3436cfSJed Brown /*@ 2901c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2911c3436cfSJed Brown 2921c3436cfSJed Brown Collective on TSAdapt 2931c3436cfSJed Brown 2941c3436cfSJed Brown Input Arguments: 2951c3436cfSJed Brown + adapt - adaptive controller context 2961c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2971c3436cfSJed Brown 2981c3436cfSJed Brown Level: intermediate 2991c3436cfSJed Brown 3001c3436cfSJed Brown .seealso: TSAdaptChoose() 3011c3436cfSJed Brown @*/ 3021c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3031c3436cfSJed Brown { 3041c3436cfSJed Brown PetscErrorCode ierr; 3051c3436cfSJed Brown 3061c3436cfSJed Brown PetscFunctionBegin; 3074782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3084782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3091c3436cfSJed Brown if (flg) { 310ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3111c3436cfSJed Brown } else { 3121c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3131c3436cfSJed Brown } 3141c3436cfSJed Brown PetscFunctionReturn(0); 3151c3436cfSJed Brown } 3161c3436cfSJed Brown 3171c3436cfSJed Brown #undef __FUNCT__ 3180873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage" 3190873808bSJed Brown /*@C 3200873808bSJed Brown TSAdaptSetCheckStage - set a callback to check convergence for a stage 3210873808bSJed Brown 3220873808bSJed Brown Logically collective on TSAdapt 3230873808bSJed Brown 3240873808bSJed Brown Input Arguments: 3250873808bSJed Brown + adapt - adaptive controller context 3260873808bSJed Brown - func - stage check function 3270873808bSJed Brown 3280873808bSJed Brown Arguments of func: 3290873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3300873808bSJed Brown 3310873808bSJed Brown + adapt - adaptive controller context 3320873808bSJed Brown . ts - time stepping context 3330873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3340873808bSJed Brown 3350873808bSJed Brown Level: advanced 3360873808bSJed Brown 3370873808bSJed Brown .seealso: TSAdaptChoose() 3380873808bSJed Brown @*/ 3390873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 3400873808bSJed Brown { 3410873808bSJed Brown 3420873808bSJed Brown PetscFunctionBegin; 3430873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3440873808bSJed Brown adapt->ops->checkstage = func; 3450873808bSJed Brown PetscFunctionReturn(0); 3460873808bSJed Brown } 3470873808bSJed Brown 3480873808bSJed Brown #undef __FUNCT__ 3491c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 3501c3436cfSJed Brown /*@ 3511c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 3521c3436cfSJed Brown 3531c3436cfSJed Brown Logically Collective 3541c3436cfSJed Brown 3551c3436cfSJed Brown Input Arguments: 356552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 3571c3436cfSJed Brown . hmin - minimum time step 3581c3436cfSJed Brown - hmax - maximum time step 3591c3436cfSJed Brown 3601c3436cfSJed Brown Options Database Keys: 3611c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 3621c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 3631c3436cfSJed Brown 3641c3436cfSJed Brown Level: intermediate 3651c3436cfSJed Brown 3661c3436cfSJed Brown .seealso: TSAdapt 3671c3436cfSJed Brown @*/ 3681c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 3691c3436cfSJed Brown { 3701c3436cfSJed Brown 3711c3436cfSJed Brown PetscFunctionBegin; 3724782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3731c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 3741c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 3751c3436cfSJed Brown PetscFunctionReturn(0); 3761c3436cfSJed Brown } 3771c3436cfSJed Brown 3781c3436cfSJed Brown #undef __FUNCT__ 37984df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 380e55864a3SBarry Smith /* 38184df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 38284df9cb4SJed Brown 38384df9cb4SJed Brown Collective on TSAdapt 38484df9cb4SJed Brown 38584df9cb4SJed Brown Input Parameter: 38684df9cb4SJed Brown . adapt - the TSAdapt context 38784df9cb4SJed Brown 38884df9cb4SJed Brown Options Database Keys: 38984df9cb4SJed Brown . -ts_adapt_type <type> - basic 39084df9cb4SJed Brown 39184df9cb4SJed Brown Level: advanced 39284df9cb4SJed Brown 39384df9cb4SJed Brown Notes: 39484df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 39584df9cb4SJed Brown 396552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType() 39784df9cb4SJed Brown 39884df9cb4SJed Brown .seealso: TSGetType() 399e55864a3SBarry Smith */ 4008c34d3f5SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptions *PetscOptionsObject,TSAdapt adapt) 40184df9cb4SJed Brown { 40284df9cb4SJed Brown PetscErrorCode ierr; 40384df9cb4SJed Brown char type[256] = TSADAPTBASIC; 4041c3436cfSJed Brown PetscBool set,flg; 40584df9cb4SJed Brown 40684df9cb4SJed Brown PetscFunctionBegin; 4074782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 40884df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 40984df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 410e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 411a264d7a6SBarry Smith ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 4128caf3d72SBarry Smith ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 41384df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 41484df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 41584df9cb4SJed Brown } 4160298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr); 4170298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr); 4180298fd71SBarry 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); 4191c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 4207619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 4217619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 4221c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 423e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 42484df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 42584df9cb4SJed Brown PetscFunctionReturn(0); 42684df9cb4SJed Brown } 42784df9cb4SJed Brown 42884df9cb4SJed Brown #undef __FUNCT__ 42984df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 43084df9cb4SJed Brown /*@ 43184df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 43284df9cb4SJed Brown 43384df9cb4SJed Brown Logically Collective 43484df9cb4SJed Brown 43584df9cb4SJed Brown Input Argument: 43684df9cb4SJed Brown . adapt - adaptive controller 43784df9cb4SJed Brown 43884df9cb4SJed Brown Level: developer 43984df9cb4SJed Brown 44084df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 44184df9cb4SJed Brown @*/ 44284df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 44384df9cb4SJed Brown { 44484df9cb4SJed Brown PetscErrorCode ierr; 44584df9cb4SJed Brown 44684df9cb4SJed Brown PetscFunctionBegin; 4474782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 44884df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 44984df9cb4SJed Brown PetscFunctionReturn(0); 45084df9cb4SJed Brown } 45184df9cb4SJed Brown 45284df9cb4SJed Brown #undef __FUNCT__ 45384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 45484df9cb4SJed Brown /*@C 45584df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 45684df9cb4SJed Brown 45784df9cb4SJed Brown Logically Collective 45884df9cb4SJed Brown 45984df9cb4SJed Brown Input Arguments: 460552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 46184df9cb4SJed Brown . name - name of the candidate scheme to add 46284df9cb4SJed Brown . order - order of the candidate scheme 46384df9cb4SJed Brown . stageorder - stage order of the candidate scheme 4648d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 46584df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 46684df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 46784df9cb4SJed Brown 46884df9cb4SJed Brown Note: 46984df9cb4SJed Brown This routine is not available in Fortran. 47084df9cb4SJed Brown 47184df9cb4SJed Brown Level: developer 47284df9cb4SJed Brown 47384df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 47484df9cb4SJed Brown @*/ 4758d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 47684df9cb4SJed Brown { 47784df9cb4SJed Brown PetscInt c; 47884df9cb4SJed Brown 47984df9cb4SJed Brown PetscFunctionBegin; 48084df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 481ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 48284df9cb4SJed Brown if (inuse) { 483ce94432eSBarry 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()"); 48484df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 48584df9cb4SJed Brown } 4861c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 4871c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 488bbd56ea5SKarl Rupp 48984df9cb4SJed Brown adapt->candidates.name[c] = name; 49084df9cb4SJed Brown adapt->candidates.order[c] = order; 49184df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 4928d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 49384df9cb4SJed Brown adapt->candidates.cost[c] = cost; 49484df9cb4SJed Brown adapt->candidates.n++; 49584df9cb4SJed Brown PetscFunctionReturn(0); 49684df9cb4SJed Brown } 49784df9cb4SJed Brown 49884df9cb4SJed Brown #undef __FUNCT__ 4998d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet" 5008d59e960SJed Brown /*@C 5018d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 5028d59e960SJed Brown 5038d59e960SJed Brown Not Collective 5048d59e960SJed Brown 5058d59e960SJed Brown Input Arguments: 5068d59e960SJed Brown . adapt - time step adaptivity context 5078d59e960SJed Brown 5088d59e960SJed Brown Output Arguments: 5098d59e960SJed Brown + n - number of candidate schemes, always at least 1 5108d59e960SJed Brown . order - the order of each candidate scheme 5118d59e960SJed Brown . stageorder - the stage order of each candidate scheme 5128d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 5138d59e960SJed Brown - cost - the relative cost of each scheme 5148d59e960SJed Brown 5158d59e960SJed Brown Level: developer 5168d59e960SJed Brown 5178d59e960SJed Brown Note: 5188d59e960SJed Brown The current scheme is always returned in the first slot 5198d59e960SJed Brown 5208d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 5218d59e960SJed Brown @*/ 5228d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 5238d59e960SJed Brown { 5248d59e960SJed Brown PetscFunctionBegin; 5258d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5268d59e960SJed Brown if (n) *n = adapt->candidates.n; 5278d59e960SJed Brown if (order) *order = adapt->candidates.order; 5288d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 5298d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 5308d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 5318d59e960SJed Brown PetscFunctionReturn(0); 5328d59e960SJed Brown } 5338d59e960SJed Brown 5348d59e960SJed Brown #undef __FUNCT__ 53584df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 53684df9cb4SJed Brown /*@C 53784df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 53884df9cb4SJed Brown 53984df9cb4SJed Brown Logically Collective 54084df9cb4SJed Brown 54184df9cb4SJed Brown Input Arguments: 54284df9cb4SJed Brown + adapt - adaptive contoller 54384df9cb4SJed Brown - h - current step size 54484df9cb4SJed Brown 54584df9cb4SJed Brown Output Arguments: 54684df9cb4SJed Brown + next_sc - scheme to use for the next step 54784df9cb4SJed Brown . next_h - step size to use for the next step 54884df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 54984df9cb4SJed Brown 5501c3436cfSJed Brown Note: 5511c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 5521c3436cfSJed Brown being retried after an initial rejection. 5531c3436cfSJed Brown 55484df9cb4SJed Brown Level: developer 55584df9cb4SJed Brown 55684df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 55784df9cb4SJed Brown @*/ 55884df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 55984df9cb4SJed Brown { 56084df9cb4SJed Brown PetscErrorCode ierr; 5610b99f514SJed Brown PetscReal wlte = -1.0; 56284df9cb4SJed Brown 56384df9cb4SJed Brown PetscFunctionBegin; 56484df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 56584df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 56684df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 56784df9cb4SJed Brown PetscValidPointer(next_h,5); 56884df9cb4SJed Brown PetscValidIntPointer(accept,6); 569ce94432eSBarry Smith if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 570ce94432eSBarry Smith if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 5710b99f514SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 57249354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 57349354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 57449354f04SShri Abhyankar PetscReal max_time = ts->max_time; 57549354f04SShri Abhyankar PetscReal next_dt = 0.0; 57649354f04SShri Abhyankar if (ts->ptime + ts->time_step + *next_h >= max_time) { 57749354f04SShri Abhyankar next_dt = max_time - (ts->ptime + ts->time_step); 5788709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 57949354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 58049354f04SShri Abhyankar } 58149354f04SShri Abhyankar } 582ce94432eSBarry Smith if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1); 58357622a8eSBarry Smith if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 5841c3436cfSJed Brown 5851c3436cfSJed Brown if (adapt->monitor) { 5861c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5870b99f514SJed Brown if (wlte < 0) { 58848c19aefSLisandro Dalcin ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr); 5890b99f514SJed Brown } else { 5906de24e2aSJed Brown ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D %s t=%-11g+%10.3e wlte=%5.3g family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)wlte,((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr); 5910b99f514SJed Brown } 5921c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5931c3436cfSJed Brown } 59484df9cb4SJed Brown PetscFunctionReturn(0); 59584df9cb4SJed Brown } 59684df9cb4SJed Brown 59784df9cb4SJed Brown #undef __FUNCT__ 59897335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage" 59997335746SJed Brown /*@ 60097335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 60197335746SJed Brown 60297335746SJed Brown Collective 60397335746SJed Brown 60497335746SJed Brown Input Arguments: 60597335746SJed Brown + adapt - adaptive controller context 60697335746SJed Brown - ts - time stepper 60797335746SJed Brown 60897335746SJed Brown Output Arguments: 60997335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 61097335746SJed Brown 61197335746SJed Brown Level: developer 61297335746SJed Brown 61397335746SJed Brown .seealso: 61497335746SJed Brown @*/ 61597335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept) 61697335746SJed Brown { 61797335746SJed Brown PetscErrorCode ierr; 61897335746SJed Brown SNES snes; 61997335746SJed Brown SNESConvergedReason snesreason; 62097335746SJed Brown 62197335746SJed Brown PetscFunctionBegin; 6224782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6234782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 6244782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 62597335746SJed Brown *accept = PETSC_TRUE; 62697335746SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 62797335746SJed Brown ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 62897335746SJed Brown if (snesreason < 0) { 62997335746SJed Brown PetscReal dt,new_dt; 63097335746SJed Brown *accept = PETSC_FALSE; 63197335746SJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 6326de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 63397335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 63448c19aefSLisandro 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); 63597335746SJed Brown if (adapt->monitor) { 63697335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 63748c19aefSLisandro Dalcin 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,dt,ts->num_snes_failures);CHKERRQ(ierr); 63897335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 63997335746SJed Brown } 64097335746SJed Brown } else { 64197335746SJed Brown new_dt = dt*adapt->scale_solve_failed; 64297335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 64397335746SJed Brown if (adapt->monitor) { 64497335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 64597335746SJed Brown 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); 64697335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 64797335746SJed Brown } 64897335746SJed Brown } 64997335746SJed Brown } 6500873808bSJed Brown if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);} 65197335746SJed Brown PetscFunctionReturn(0); 65297335746SJed Brown } 65397335746SJed Brown 6540873808bSJed Brown 6550873808bSJed Brown 65697335746SJed Brown #undef __FUNCT__ 65784df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 65884df9cb4SJed Brown /*@ 65984df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 66084df9cb4SJed Brown 66184df9cb4SJed Brown Collective on MPI_Comm 66284df9cb4SJed Brown 66384df9cb4SJed Brown Input Parameter: 66484df9cb4SJed Brown . comm - The communicator 66584df9cb4SJed Brown 66684df9cb4SJed Brown Output Parameter: 66784df9cb4SJed Brown . adapt - new TSAdapt object 66884df9cb4SJed Brown 66984df9cb4SJed Brown Level: developer 67084df9cb4SJed Brown 67184df9cb4SJed Brown Notes: 67284df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 67384df9cb4SJed Brown 67484df9cb4SJed Brown .keywords: TSAdapt, create 675552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 67684df9cb4SJed Brown @*/ 67784df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 67884df9cb4SJed Brown { 67984df9cb4SJed Brown PetscErrorCode ierr; 68084df9cb4SJed Brown TSAdapt adapt; 68184df9cb4SJed Brown 68284df9cb4SJed Brown PetscFunctionBegin; 6833b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 6843b3bcf4cSLisandro Dalcin *inadapt = NULL; 6853b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 6863b3bcf4cSLisandro Dalcin 68767c2884eSBarry Smith ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 6881c3436cfSJed Brown 6891c3436cfSJed Brown adapt->dt_min = 1e-20; 6901c3436cfSJed Brown adapt->dt_max = 1e50; 69197335746SJed Brown adapt->scale_solve_failed = 0.25; 6927619abb3SShri adapt->wnormtype = NORM_2; 6931c3436cfSJed Brown 69484df9cb4SJed Brown *inadapt = adapt; 69584df9cb4SJed Brown PetscFunctionReturn(0); 69684df9cb4SJed Brown } 697