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; 70bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 71bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 72bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 7384df9cb4SJed Brown PetscFunctionReturn(0); 7484df9cb4SJed Brown } 7584df9cb4SJed Brown 7684df9cb4SJed Brown #undef __FUNCT__ 7784df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage" 7884df9cb4SJed Brown /*@C 7984df9cb4SJed Brown TSFinalizePackage - This function destroys everything in the TS package. It is 8084df9cb4SJed Brown called from PetscFinalize(). 8184df9cb4SJed Brown 8284df9cb4SJed Brown Level: developer 8384df9cb4SJed Brown 8484df9cb4SJed Brown .keywords: Petsc, destroy, package 8584df9cb4SJed Brown .seealso: PetscFinalize() 8684df9cb4SJed Brown @*/ 8784df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8884df9cb4SJed Brown { 8937e93019SBarry Smith PetscErrorCode ierr; 9037e93019SBarry Smith 9184df9cb4SJed Brown PetscFunctionBegin; 9237e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9384df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9484df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9584df9cb4SJed Brown PetscFunctionReturn(0); 9684df9cb4SJed Brown } 9784df9cb4SJed Brown 9884df9cb4SJed Brown #undef __FUNCT__ 9984df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage" 10084df9cb4SJed Brown /*@C 10184df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 10284df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 10384df9cb4SJed Brown TSCreate_GL() when using static libraries. 10484df9cb4SJed Brown 10584df9cb4SJed Brown Level: developer 10684df9cb4SJed Brown 10784df9cb4SJed Brown .keywords: TSAdapt, initialize, package 10884df9cb4SJed Brown .seealso: PetscInitialize() 10984df9cb4SJed Brown @*/ 110607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11184df9cb4SJed Brown { 11284df9cb4SJed Brown PetscErrorCode ierr; 11384df9cb4SJed Brown 11484df9cb4SJed Brown PetscFunctionBegin; 11584df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11684df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11784df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 118607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 11984df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 12084df9cb4SJed Brown PetscFunctionReturn(0); 12184df9cb4SJed Brown } 12284df9cb4SJed Brown 12384df9cb4SJed Brown #undef __FUNCT__ 12484df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 12519fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 12684df9cb4SJed Brown { 12784df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 12884df9cb4SJed Brown 12984df9cb4SJed Brown PetscFunctionBegin; 1301c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 13184df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 13284df9cb4SJed Brown if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 133*4cefc2ffSBarry Smith /* Reinitialize function pointers in TSAdaptOps structure */ 134*4cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 13584df9cb4SJed Brown ierr = (*r)(adapt);CHKERRQ(ierr); 13684df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 13784df9cb4SJed Brown PetscFunctionReturn(0); 13884df9cb4SJed Brown } 13984df9cb4SJed Brown 14084df9cb4SJed Brown #undef __FUNCT__ 14184df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 14284df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 14384df9cb4SJed Brown { 14484df9cb4SJed Brown PetscErrorCode ierr; 14584df9cb4SJed Brown 14684df9cb4SJed Brown PetscFunctionBegin; 14784df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 14884df9cb4SJed Brown PetscFunctionReturn(0); 14984df9cb4SJed Brown } 15084df9cb4SJed Brown 15184df9cb4SJed Brown #undef __FUNCT__ 152ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad" 153ad6bc421SBarry Smith /*@C 154ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 155ad6bc421SBarry Smith 156ad6bc421SBarry Smith Collective on PetscViewer 157ad6bc421SBarry Smith 158ad6bc421SBarry Smith Input Parameters: 159ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 160ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 161ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 162ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 163ad6bc421SBarry Smith 164ad6bc421SBarry Smith Level: intermediate 165ad6bc421SBarry Smith 166ad6bc421SBarry Smith Notes: 167ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 168ad6bc421SBarry Smith 169ad6bc421SBarry Smith Notes for advanced users: 170ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 171ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 172ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 173ad6bc421SBarry Smith format is 174ad6bc421SBarry Smith .vb 175ad6bc421SBarry Smith has not yet been determined 176ad6bc421SBarry Smith .ve 177ad6bc421SBarry Smith 178ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 179ad6bc421SBarry Smith @*/ 180ad6bc421SBarry Smith PetscErrorCode TSAdaptLoad(TSAdapt tsadapt, PetscViewer viewer) 181ad6bc421SBarry Smith { 182ad6bc421SBarry Smith PetscErrorCode ierr; 183ad6bc421SBarry Smith PetscBool isbinary; 184ad6bc421SBarry Smith char type[256]; 185ad6bc421SBarry Smith 186ad6bc421SBarry Smith PetscFunctionBegin; 187ad6bc421SBarry Smith PetscValidHeaderSpecific(tsadapt,TSADAPT_CLASSID,1); 188ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 189ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 190ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 191ad6bc421SBarry Smith 192ad6bc421SBarry Smith ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 193ad6bc421SBarry Smith ierr = TSAdaptSetType(tsadapt, type);CHKERRQ(ierr); 194ad6bc421SBarry Smith if (tsadapt->ops->load) { 195ad6bc421SBarry Smith ierr = (*tsadapt->ops->load)(tsadapt,viewer);CHKERRQ(ierr); 196ad6bc421SBarry Smith } 197ad6bc421SBarry Smith PetscFunctionReturn(0); 198ad6bc421SBarry Smith } 199ad6bc421SBarry Smith 200ad6bc421SBarry Smith #undef __FUNCT__ 20184df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 20284df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 20384df9cb4SJed Brown { 20484df9cb4SJed Brown PetscErrorCode ierr; 205ad6bc421SBarry Smith PetscBool iascii,isbinary; 20684df9cb4SJed Brown 20784df9cb4SJed Brown PetscFunctionBegin; 208251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 209ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 21084df9cb4SJed Brown if (iascii) { 211dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 21284df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 21384df9cb4SJed Brown if (adapt->ops->view) { 21484df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 21584df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 21684df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 21784df9cb4SJed Brown } 218ad6bc421SBarry Smith } else if (isbinary) { 219ad6bc421SBarry Smith char type[256]; 220ad6bc421SBarry Smith 221ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 222ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 223ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 224bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 225f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 226f2c2a1b9SBarry Smith } 22784df9cb4SJed Brown PetscFunctionReturn(0); 22884df9cb4SJed Brown } 22984df9cb4SJed Brown 23084df9cb4SJed Brown #undef __FUNCT__ 23184df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 23284df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 23384df9cb4SJed Brown { 23484df9cb4SJed Brown PetscErrorCode ierr; 23584df9cb4SJed Brown 23684df9cb4SJed Brown PetscFunctionBegin; 23784df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 23884df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 23984df9cb4SJed Brown if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 24084df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 2411c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 24284df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 24384df9cb4SJed Brown PetscFunctionReturn(0); 24484df9cb4SJed Brown } 24584df9cb4SJed Brown 24684df9cb4SJed Brown #undef __FUNCT__ 2471c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 2481c3436cfSJed Brown /*@ 2491c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2501c3436cfSJed Brown 2511c3436cfSJed Brown Collective on TSAdapt 2521c3436cfSJed Brown 2531c3436cfSJed Brown Input Arguments: 2541c3436cfSJed Brown + adapt - adaptive controller context 2551c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2561c3436cfSJed Brown 2571c3436cfSJed Brown Level: intermediate 2581c3436cfSJed Brown 2591c3436cfSJed Brown .seealso: TSAdaptChoose() 2601c3436cfSJed Brown @*/ 2611c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 2621c3436cfSJed Brown { 2631c3436cfSJed Brown PetscErrorCode ierr; 2641c3436cfSJed Brown 2651c3436cfSJed Brown PetscFunctionBegin; 2661c3436cfSJed Brown if (flg) { 267ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 2681c3436cfSJed Brown } else { 2691c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 2701c3436cfSJed Brown } 2711c3436cfSJed Brown PetscFunctionReturn(0); 2721c3436cfSJed Brown } 2731c3436cfSJed Brown 2741c3436cfSJed Brown #undef __FUNCT__ 2750873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage" 2760873808bSJed Brown /*@C 2770873808bSJed Brown TSAdaptSetCheckStage - set a callback to check convergence for a stage 2780873808bSJed Brown 2790873808bSJed Brown Logically collective on TSAdapt 2800873808bSJed Brown 2810873808bSJed Brown Input Arguments: 2820873808bSJed Brown + adapt - adaptive controller context 2830873808bSJed Brown - func - stage check function 2840873808bSJed Brown 2850873808bSJed Brown Arguments of func: 2860873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 2870873808bSJed Brown 2880873808bSJed Brown + adapt - adaptive controller context 2890873808bSJed Brown . ts - time stepping context 2900873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 2910873808bSJed Brown 2920873808bSJed Brown Level: advanced 2930873808bSJed Brown 2940873808bSJed Brown .seealso: TSAdaptChoose() 2950873808bSJed Brown @*/ 2960873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 2970873808bSJed Brown { 2980873808bSJed Brown 2990873808bSJed Brown PetscFunctionBegin; 3000873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3010873808bSJed Brown adapt->ops->checkstage = func; 3020873808bSJed Brown PetscFunctionReturn(0); 3030873808bSJed Brown } 3040873808bSJed Brown 3050873808bSJed Brown #undef __FUNCT__ 3061c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 3071c3436cfSJed Brown /*@ 3081c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 3091c3436cfSJed Brown 3101c3436cfSJed Brown Logically Collective 3111c3436cfSJed Brown 3121c3436cfSJed Brown Input Arguments: 313552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 3141c3436cfSJed Brown . hmin - minimum time step 3151c3436cfSJed Brown - hmax - maximum time step 3161c3436cfSJed Brown 3171c3436cfSJed Brown Options Database Keys: 3181c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 3191c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 3201c3436cfSJed Brown 3211c3436cfSJed Brown Level: intermediate 3221c3436cfSJed Brown 3231c3436cfSJed Brown .seealso: TSAdapt 3241c3436cfSJed Brown @*/ 3251c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 3261c3436cfSJed Brown { 3271c3436cfSJed Brown 3281c3436cfSJed Brown PetscFunctionBegin; 3291c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 3301c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 3311c3436cfSJed Brown PetscFunctionReturn(0); 3321c3436cfSJed Brown } 3331c3436cfSJed Brown 3341c3436cfSJed Brown #undef __FUNCT__ 33584df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 33684df9cb4SJed Brown /*@ 33784df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 33884df9cb4SJed Brown 33984df9cb4SJed Brown Collective on TSAdapt 34084df9cb4SJed Brown 34184df9cb4SJed Brown Input Parameter: 34284df9cb4SJed Brown . adapt - the TSAdapt context 34384df9cb4SJed Brown 34484df9cb4SJed Brown Options Database Keys: 34584df9cb4SJed Brown . -ts_adapt_type <type> - basic 34684df9cb4SJed Brown 34784df9cb4SJed Brown Level: advanced 34884df9cb4SJed Brown 34984df9cb4SJed Brown Notes: 35084df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 35184df9cb4SJed Brown 352552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType() 35384df9cb4SJed Brown 35484df9cb4SJed Brown .seealso: TSGetType() 35584df9cb4SJed Brown @*/ 35684df9cb4SJed Brown PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt) 35784df9cb4SJed Brown { 35884df9cb4SJed Brown PetscErrorCode ierr; 35984df9cb4SJed Brown char type[256] = TSADAPTBASIC; 3601c3436cfSJed Brown PetscBool set,flg; 36184df9cb4SJed Brown 36284df9cb4SJed Brown PetscFunctionBegin; 36384df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 36484df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 36584df9cb4SJed Brown ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr); 366a264d7a6SBarry Smith ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 3678caf3d72SBarry Smith ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 36884df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 36984df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 37084df9cb4SJed Brown } 37184df9cb4SJed Brown if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);} 3720298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr); 3730298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr); 3740298fd71SBarry 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); 3751c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 3761c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 37784df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 37884df9cb4SJed Brown PetscFunctionReturn(0); 37984df9cb4SJed Brown } 38084df9cb4SJed Brown 38184df9cb4SJed Brown #undef __FUNCT__ 38284df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 38384df9cb4SJed Brown /*@ 38484df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 38584df9cb4SJed Brown 38684df9cb4SJed Brown Logically Collective 38784df9cb4SJed Brown 38884df9cb4SJed Brown Input Argument: 38984df9cb4SJed Brown . adapt - adaptive controller 39084df9cb4SJed Brown 39184df9cb4SJed Brown Level: developer 39284df9cb4SJed Brown 39384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 39484df9cb4SJed Brown @*/ 39584df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 39684df9cb4SJed Brown { 39784df9cb4SJed Brown PetscErrorCode ierr; 39884df9cb4SJed Brown 39984df9cb4SJed Brown PetscFunctionBegin; 40084df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 40184df9cb4SJed Brown PetscFunctionReturn(0); 40284df9cb4SJed Brown } 40384df9cb4SJed Brown 40484df9cb4SJed Brown #undef __FUNCT__ 40584df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 40684df9cb4SJed Brown /*@C 40784df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 40884df9cb4SJed Brown 40984df9cb4SJed Brown Logically Collective 41084df9cb4SJed Brown 41184df9cb4SJed Brown Input Arguments: 412552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 41384df9cb4SJed Brown . name - name of the candidate scheme to add 41484df9cb4SJed Brown . order - order of the candidate scheme 41584df9cb4SJed Brown . stageorder - stage order of the candidate scheme 4168d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 41784df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 41884df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 41984df9cb4SJed Brown 42084df9cb4SJed Brown Note: 42184df9cb4SJed Brown This routine is not available in Fortran. 42284df9cb4SJed Brown 42384df9cb4SJed Brown Level: developer 42484df9cb4SJed Brown 42584df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 42684df9cb4SJed Brown @*/ 4278d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 42884df9cb4SJed Brown { 42984df9cb4SJed Brown PetscInt c; 43084df9cb4SJed Brown 43184df9cb4SJed Brown PetscFunctionBegin; 43284df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 433ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 43484df9cb4SJed Brown if (inuse) { 435ce94432eSBarry 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()"); 43684df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 43784df9cb4SJed Brown } 4381c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 4391c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 440bbd56ea5SKarl Rupp 44184df9cb4SJed Brown adapt->candidates.name[c] = name; 44284df9cb4SJed Brown adapt->candidates.order[c] = order; 44384df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 4448d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 44584df9cb4SJed Brown adapt->candidates.cost[c] = cost; 44684df9cb4SJed Brown adapt->candidates.n++; 44784df9cb4SJed Brown PetscFunctionReturn(0); 44884df9cb4SJed Brown } 44984df9cb4SJed Brown 45084df9cb4SJed Brown #undef __FUNCT__ 4518d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet" 4528d59e960SJed Brown /*@C 4538d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 4548d59e960SJed Brown 4558d59e960SJed Brown Not Collective 4568d59e960SJed Brown 4578d59e960SJed Brown Input Arguments: 4588d59e960SJed Brown . adapt - time step adaptivity context 4598d59e960SJed Brown 4608d59e960SJed Brown Output Arguments: 4618d59e960SJed Brown + n - number of candidate schemes, always at least 1 4628d59e960SJed Brown . order - the order of each candidate scheme 4638d59e960SJed Brown . stageorder - the stage order of each candidate scheme 4648d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 4658d59e960SJed Brown - cost - the relative cost of each scheme 4668d59e960SJed Brown 4678d59e960SJed Brown Level: developer 4688d59e960SJed Brown 4698d59e960SJed Brown Note: 4708d59e960SJed Brown The current scheme is always returned in the first slot 4718d59e960SJed Brown 4728d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 4738d59e960SJed Brown @*/ 4748d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 4758d59e960SJed Brown { 4768d59e960SJed Brown PetscFunctionBegin; 4778d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4788d59e960SJed Brown if (n) *n = adapt->candidates.n; 4798d59e960SJed Brown if (order) *order = adapt->candidates.order; 4808d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 4818d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 4828d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 4838d59e960SJed Brown PetscFunctionReturn(0); 4848d59e960SJed Brown } 4858d59e960SJed Brown 4868d59e960SJed Brown #undef __FUNCT__ 48784df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 48884df9cb4SJed Brown /*@C 48984df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 49084df9cb4SJed Brown 49184df9cb4SJed Brown Logically Collective 49284df9cb4SJed Brown 49384df9cb4SJed Brown Input Arguments: 49484df9cb4SJed Brown + adapt - adaptive contoller 49584df9cb4SJed Brown - h - current step size 49684df9cb4SJed Brown 49784df9cb4SJed Brown Output Arguments: 49884df9cb4SJed Brown + next_sc - scheme to use for the next step 49984df9cb4SJed Brown . next_h - step size to use for the next step 50084df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 50184df9cb4SJed Brown 5021c3436cfSJed Brown Note: 5031c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 5041c3436cfSJed Brown being retried after an initial rejection. 5051c3436cfSJed Brown 50684df9cb4SJed Brown Level: developer 50784df9cb4SJed Brown 50884df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 50984df9cb4SJed Brown @*/ 51084df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 51184df9cb4SJed Brown { 51284df9cb4SJed Brown PetscErrorCode ierr; 5130b99f514SJed Brown PetscReal wlte = -1.0; 51484df9cb4SJed Brown 51584df9cb4SJed Brown PetscFunctionBegin; 51684df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 51784df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 51884df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 51984df9cb4SJed Brown PetscValidPointer(next_h,5); 52084df9cb4SJed Brown PetscValidIntPointer(accept,6); 521ce94432eSBarry Smith if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 522ce94432eSBarry 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); 5230b99f514SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 52449354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 52549354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 52649354f04SShri Abhyankar PetscReal max_time = ts->max_time; 52749354f04SShri Abhyankar PetscReal next_dt = 0.0; 52849354f04SShri Abhyankar if (ts->ptime + ts->time_step + *next_h >= max_time) { 52949354f04SShri Abhyankar next_dt = max_time - (ts->ptime + ts->time_step); 5308709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 53149354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 53249354f04SShri Abhyankar } 53349354f04SShri Abhyankar } 534ce94432eSBarry 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); 53557622a8eSBarry Smith if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 5361c3436cfSJed Brown 5371c3436cfSJed Brown if (adapt->monitor) { 5381c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5390b99f514SJed Brown if (wlte < 0) { 5406de24e2aSJed Brown ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10g\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); 5410b99f514SJed Brown } else { 5426de24e2aSJed 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); 5430b99f514SJed Brown } 5441c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5451c3436cfSJed Brown } 54684df9cb4SJed Brown PetscFunctionReturn(0); 54784df9cb4SJed Brown } 54884df9cb4SJed Brown 54984df9cb4SJed Brown #undef __FUNCT__ 55097335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage" 55197335746SJed Brown /*@ 55297335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 55397335746SJed Brown 55497335746SJed Brown Collective 55597335746SJed Brown 55697335746SJed Brown Input Arguments: 55797335746SJed Brown + adapt - adaptive controller context 55897335746SJed Brown - ts - time stepper 55997335746SJed Brown 56097335746SJed Brown Output Arguments: 56197335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 56297335746SJed Brown 56397335746SJed Brown Level: developer 56497335746SJed Brown 56597335746SJed Brown .seealso: 56697335746SJed Brown @*/ 56797335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept) 56897335746SJed Brown { 56997335746SJed Brown PetscErrorCode ierr; 57097335746SJed Brown SNES snes; 57197335746SJed Brown SNESConvergedReason snesreason; 57297335746SJed Brown 57397335746SJed Brown PetscFunctionBegin; 57497335746SJed Brown *accept = PETSC_TRUE; 57597335746SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 57697335746SJed Brown ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 57797335746SJed Brown if (snesreason < 0) { 57897335746SJed Brown PetscReal dt,new_dt; 57997335746SJed Brown *accept = PETSC_FALSE; 58097335746SJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 5816de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 58297335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 58397335746SJed Brown ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 58497335746SJed Brown if (adapt->monitor) { 58597335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 58697335746SJed Brown ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e, %D failures exceeds current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,dt,ts->num_snes_failures);CHKERRQ(ierr); 58797335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 58897335746SJed Brown } 58997335746SJed Brown } else { 59097335746SJed Brown new_dt = dt*adapt->scale_solve_failed; 59197335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 59297335746SJed Brown if (adapt->monitor) { 59397335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 59497335746SJed 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); 59597335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 59697335746SJed Brown } 59797335746SJed Brown } 59897335746SJed Brown } 5990873808bSJed Brown if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);} 60097335746SJed Brown PetscFunctionReturn(0); 60197335746SJed Brown } 60297335746SJed Brown 6030873808bSJed Brown 6040873808bSJed Brown 60597335746SJed Brown #undef __FUNCT__ 60684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 60784df9cb4SJed Brown /*@ 60884df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 60984df9cb4SJed Brown 61084df9cb4SJed Brown Collective on MPI_Comm 61184df9cb4SJed Brown 61284df9cb4SJed Brown Input Parameter: 61384df9cb4SJed Brown . comm - The communicator 61484df9cb4SJed Brown 61584df9cb4SJed Brown Output Parameter: 61684df9cb4SJed Brown . adapt - new TSAdapt object 61784df9cb4SJed Brown 61884df9cb4SJed Brown Level: developer 61984df9cb4SJed Brown 62084df9cb4SJed Brown Notes: 62184df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 62284df9cb4SJed Brown 62384df9cb4SJed Brown .keywords: TSAdapt, create 624552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 62584df9cb4SJed Brown @*/ 62684df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 62784df9cb4SJed Brown { 62884df9cb4SJed Brown PetscErrorCode ierr; 62984df9cb4SJed Brown TSAdapt adapt; 63084df9cb4SJed Brown 63184df9cb4SJed Brown PetscFunctionBegin; 63284df9cb4SJed Brown *inadapt = 0; 63367c2884eSBarry Smith ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 6341c3436cfSJed Brown 6351c3436cfSJed Brown adapt->dt_min = 1e-20; 6361c3436cfSJed Brown adapt->dt_max = 1e50; 63797335746SJed Brown adapt->scale_solve_failed = 0.25; 6381c3436cfSJed Brown 63984df9cb4SJed Brown *inadapt = adapt; 64084df9cb4SJed Brown PetscFunctionReturn(0); 64184df9cb4SJed Brown } 642