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 { 89*37e93019SBarry Smith PetscErrorCode ierr; 90*37e93019SBarry Smith 9184df9cb4SJed Brown PetscFunctionBegin; 92*37e93019SBarry 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__ "TSAdaptRegisterDestroy" 12584df9cb4SJed Brown /*@C 126607a6623SBarry Smith TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister() 12784df9cb4SJed Brown 12884df9cb4SJed Brown Not Collective 12984df9cb4SJed Brown 13084df9cb4SJed Brown Level: advanced 13184df9cb4SJed Brown 13284df9cb4SJed Brown .keywords: TSAdapt, register, destroy 133607a6623SBarry Smith .seealso: TSAdaptRegister(), TSAdaptRegisterAll() 13484df9cb4SJed Brown @*/ 13584df9cb4SJed Brown PetscErrorCode TSAdaptRegisterDestroy(void) 13684df9cb4SJed Brown { 13784df9cb4SJed Brown PetscErrorCode ierr; 13884df9cb4SJed Brown 13984df9cb4SJed Brown PetscFunctionBegin; 140140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 14184df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 14284df9cb4SJed Brown PetscFunctionReturn(0); 14384df9cb4SJed Brown } 14484df9cb4SJed Brown 14584df9cb4SJed Brown 14684df9cb4SJed Brown #undef __FUNCT__ 14784df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 14819fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 14984df9cb4SJed Brown { 15084df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 15184df9cb4SJed Brown 15284df9cb4SJed Brown PetscFunctionBegin; 1531c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 15484df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 15584df9cb4SJed Brown if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 15684df9cb4SJed Brown ierr = (*r)(adapt);CHKERRQ(ierr); 15784df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15884df9cb4SJed Brown PetscFunctionReturn(0); 15984df9cb4SJed Brown } 16084df9cb4SJed Brown 16184df9cb4SJed Brown #undef __FUNCT__ 16284df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 16384df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 16484df9cb4SJed Brown { 16584df9cb4SJed Brown PetscErrorCode ierr; 16684df9cb4SJed Brown 16784df9cb4SJed Brown PetscFunctionBegin; 16884df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 16984df9cb4SJed Brown PetscFunctionReturn(0); 17084df9cb4SJed Brown } 17184df9cb4SJed Brown 17284df9cb4SJed Brown #undef __FUNCT__ 173ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad" 174ad6bc421SBarry Smith /*@C 175ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 176ad6bc421SBarry Smith 177ad6bc421SBarry Smith Collective on PetscViewer 178ad6bc421SBarry Smith 179ad6bc421SBarry Smith Input Parameters: 180ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 181ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 182ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 183ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 184ad6bc421SBarry Smith 185ad6bc421SBarry Smith Level: intermediate 186ad6bc421SBarry Smith 187ad6bc421SBarry Smith Notes: 188ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 189ad6bc421SBarry Smith 190ad6bc421SBarry Smith Notes for advanced users: 191ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 192ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 193ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 194ad6bc421SBarry Smith format is 195ad6bc421SBarry Smith .vb 196ad6bc421SBarry Smith has not yet been determined 197ad6bc421SBarry Smith .ve 198ad6bc421SBarry Smith 199ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 200ad6bc421SBarry Smith @*/ 201ad6bc421SBarry Smith PetscErrorCode TSAdaptLoad(TSAdapt tsadapt, PetscViewer viewer) 202ad6bc421SBarry Smith { 203ad6bc421SBarry Smith PetscErrorCode ierr; 204ad6bc421SBarry Smith PetscBool isbinary; 205ad6bc421SBarry Smith char type[256]; 206ad6bc421SBarry Smith 207ad6bc421SBarry Smith PetscFunctionBegin; 208ad6bc421SBarry Smith PetscValidHeaderSpecific(tsadapt,TSADAPT_CLASSID,1); 209ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 210ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 211ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 212ad6bc421SBarry Smith 213ad6bc421SBarry Smith ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 214ad6bc421SBarry Smith ierr = TSAdaptSetType(tsadapt, type);CHKERRQ(ierr); 215ad6bc421SBarry Smith if (tsadapt->ops->load) { 216ad6bc421SBarry Smith ierr = (*tsadapt->ops->load)(tsadapt,viewer);CHKERRQ(ierr); 217ad6bc421SBarry Smith } 218ad6bc421SBarry Smith PetscFunctionReturn(0); 219ad6bc421SBarry Smith } 220ad6bc421SBarry Smith 221ad6bc421SBarry Smith #undef __FUNCT__ 22284df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 22384df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 22484df9cb4SJed Brown { 22584df9cb4SJed Brown PetscErrorCode ierr; 226ad6bc421SBarry Smith PetscBool iascii,isbinary; 22784df9cb4SJed Brown 22884df9cb4SJed Brown PetscFunctionBegin; 229251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 230ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 23184df9cb4SJed Brown if (iascii) { 23284df9cb4SJed Brown ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr); 23384df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 23484df9cb4SJed Brown if (adapt->ops->view) { 23584df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 23684df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 23784df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 23884df9cb4SJed Brown } 239ad6bc421SBarry Smith } else if (isbinary) { 240ad6bc421SBarry Smith char type[256]; 241ad6bc421SBarry Smith 242ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 243ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 244ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 245bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 246f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 247f2c2a1b9SBarry Smith } 24884df9cb4SJed Brown PetscFunctionReturn(0); 24984df9cb4SJed Brown } 25084df9cb4SJed Brown 25184df9cb4SJed Brown #undef __FUNCT__ 25284df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 25384df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 25484df9cb4SJed Brown { 25584df9cb4SJed Brown PetscErrorCode ierr; 25684df9cb4SJed Brown 25784df9cb4SJed Brown PetscFunctionBegin; 25884df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 25984df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 26084df9cb4SJed Brown if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 26184df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 2621c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 26384df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 26484df9cb4SJed Brown PetscFunctionReturn(0); 26584df9cb4SJed Brown } 26684df9cb4SJed Brown 26784df9cb4SJed Brown #undef __FUNCT__ 2681c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 2691c3436cfSJed Brown /*@ 2701c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2711c3436cfSJed Brown 2721c3436cfSJed Brown Collective on TSAdapt 2731c3436cfSJed Brown 2741c3436cfSJed Brown Input Arguments: 2751c3436cfSJed Brown + adapt - adaptive controller context 2761c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2771c3436cfSJed Brown 2781c3436cfSJed Brown Level: intermediate 2791c3436cfSJed Brown 2801c3436cfSJed Brown .seealso: TSAdaptChoose() 2811c3436cfSJed Brown @*/ 2821c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 2831c3436cfSJed Brown { 2841c3436cfSJed Brown PetscErrorCode ierr; 2851c3436cfSJed Brown 2861c3436cfSJed Brown PetscFunctionBegin; 2871c3436cfSJed Brown if (flg) { 288ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 2891c3436cfSJed Brown } else { 2901c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 2911c3436cfSJed Brown } 2921c3436cfSJed Brown PetscFunctionReturn(0); 2931c3436cfSJed Brown } 2941c3436cfSJed Brown 2951c3436cfSJed Brown #undef __FUNCT__ 2960873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage" 2970873808bSJed Brown /*@C 2980873808bSJed Brown TSAdaptSetCheckStage - set a callback to check convergence for a stage 2990873808bSJed Brown 3000873808bSJed Brown Logically collective on TSAdapt 3010873808bSJed Brown 3020873808bSJed Brown Input Arguments: 3030873808bSJed Brown + adapt - adaptive controller context 3040873808bSJed Brown - func - stage check function 3050873808bSJed Brown 3060873808bSJed Brown Arguments of func: 3070873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3080873808bSJed Brown 3090873808bSJed Brown + adapt - adaptive controller context 3100873808bSJed Brown . ts - time stepping context 3110873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3120873808bSJed Brown 3130873808bSJed Brown Level: advanced 3140873808bSJed Brown 3150873808bSJed Brown .seealso: TSAdaptChoose() 3160873808bSJed Brown @*/ 3170873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 3180873808bSJed Brown { 3190873808bSJed Brown 3200873808bSJed Brown PetscFunctionBegin; 3210873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3220873808bSJed Brown adapt->ops->checkstage = func; 3230873808bSJed Brown PetscFunctionReturn(0); 3240873808bSJed Brown } 3250873808bSJed Brown 3260873808bSJed Brown #undef __FUNCT__ 3271c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 3281c3436cfSJed Brown /*@ 3291c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 3301c3436cfSJed Brown 3311c3436cfSJed Brown Logically Collective 3321c3436cfSJed Brown 3331c3436cfSJed Brown Input Arguments: 334ad6bc421SBarry Smith + adapt - time step adaptivity context, usually gotten with TSGetTSAdapt() 3351c3436cfSJed Brown . hmin - minimum time step 3361c3436cfSJed Brown - hmax - maximum time step 3371c3436cfSJed Brown 3381c3436cfSJed Brown Options Database Keys: 3391c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 3401c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 3411c3436cfSJed Brown 3421c3436cfSJed Brown Level: intermediate 3431c3436cfSJed Brown 3441c3436cfSJed Brown .seealso: TSAdapt 3451c3436cfSJed Brown @*/ 3461c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 3471c3436cfSJed Brown { 3481c3436cfSJed Brown 3491c3436cfSJed Brown PetscFunctionBegin; 3501c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 3511c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 3521c3436cfSJed Brown PetscFunctionReturn(0); 3531c3436cfSJed Brown } 3541c3436cfSJed Brown 3551c3436cfSJed Brown #undef __FUNCT__ 35684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 35784df9cb4SJed Brown /*@ 35884df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 35984df9cb4SJed Brown 36084df9cb4SJed Brown Collective on TSAdapt 36184df9cb4SJed Brown 36284df9cb4SJed Brown Input Parameter: 36384df9cb4SJed Brown . adapt - the TSAdapt context 36484df9cb4SJed Brown 36584df9cb4SJed Brown Options Database Keys: 36684df9cb4SJed Brown . -ts_adapt_type <type> - basic 36784df9cb4SJed Brown 36884df9cb4SJed Brown Level: advanced 36984df9cb4SJed Brown 37084df9cb4SJed Brown Notes: 37184df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 37284df9cb4SJed Brown 373ad6bc421SBarry Smith .keywords: TS, TSGetTSAdapt(), TSAdaptSetType() 37484df9cb4SJed Brown 37584df9cb4SJed Brown .seealso: TSGetType() 37684df9cb4SJed Brown @*/ 37784df9cb4SJed Brown PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt) 37884df9cb4SJed Brown { 37984df9cb4SJed Brown PetscErrorCode ierr; 38084df9cb4SJed Brown char type[256] = TSADAPTBASIC; 3811c3436cfSJed Brown PetscBool set,flg; 38284df9cb4SJed Brown 38384df9cb4SJed Brown PetscFunctionBegin; 38484df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 38584df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 38684df9cb4SJed Brown ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr); 38784df9cb4SJed Brown ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 3888caf3d72SBarry Smith ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 38984df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 39084df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 39184df9cb4SJed Brown } 39284df9cb4SJed Brown if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);} 3930298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr); 3940298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr); 3950298fd71SBarry 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); 3961c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 3971c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 39884df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 39984df9cb4SJed Brown PetscFunctionReturn(0); 40084df9cb4SJed Brown } 40184df9cb4SJed Brown 40284df9cb4SJed Brown #undef __FUNCT__ 40384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 40484df9cb4SJed Brown /*@ 40584df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 40684df9cb4SJed Brown 40784df9cb4SJed Brown Logically Collective 40884df9cb4SJed Brown 40984df9cb4SJed Brown Input Argument: 41084df9cb4SJed Brown . adapt - adaptive controller 41184df9cb4SJed Brown 41284df9cb4SJed Brown Level: developer 41384df9cb4SJed Brown 41484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 41584df9cb4SJed Brown @*/ 41684df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 41784df9cb4SJed Brown { 41884df9cb4SJed Brown PetscErrorCode ierr; 41984df9cb4SJed Brown 42084df9cb4SJed Brown PetscFunctionBegin; 42184df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 42284df9cb4SJed Brown PetscFunctionReturn(0); 42384df9cb4SJed Brown } 42484df9cb4SJed Brown 42584df9cb4SJed Brown #undef __FUNCT__ 42684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 42784df9cb4SJed Brown /*@C 42884df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 42984df9cb4SJed Brown 43084df9cb4SJed Brown Logically Collective 43184df9cb4SJed Brown 43284df9cb4SJed Brown Input Arguments: 433ad6bc421SBarry Smith + adapt - time step adaptivity context, obtained with TSGetTSAdapt() or TSAdaptCreate() 43484df9cb4SJed Brown . name - name of the candidate scheme to add 43584df9cb4SJed Brown . order - order of the candidate scheme 43684df9cb4SJed Brown . stageorder - stage order of the candidate scheme 4378d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 43884df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 43984df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 44084df9cb4SJed Brown 44184df9cb4SJed Brown Note: 44284df9cb4SJed Brown This routine is not available in Fortran. 44384df9cb4SJed Brown 44484df9cb4SJed Brown Level: developer 44584df9cb4SJed Brown 44684df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 44784df9cb4SJed Brown @*/ 4488d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 44984df9cb4SJed Brown { 45084df9cb4SJed Brown PetscInt c; 45184df9cb4SJed Brown 45284df9cb4SJed Brown PetscFunctionBegin; 45384df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 454ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 45584df9cb4SJed Brown if (inuse) { 456ce94432eSBarry 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()"); 45784df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 45884df9cb4SJed Brown } 4591c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 4601c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 461bbd56ea5SKarl Rupp 46284df9cb4SJed Brown adapt->candidates.name[c] = name; 46384df9cb4SJed Brown adapt->candidates.order[c] = order; 46484df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 4658d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 46684df9cb4SJed Brown adapt->candidates.cost[c] = cost; 46784df9cb4SJed Brown adapt->candidates.n++; 46884df9cb4SJed Brown PetscFunctionReturn(0); 46984df9cb4SJed Brown } 47084df9cb4SJed Brown 47184df9cb4SJed Brown #undef __FUNCT__ 4728d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet" 4738d59e960SJed Brown /*@C 4748d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 4758d59e960SJed Brown 4768d59e960SJed Brown Not Collective 4778d59e960SJed Brown 4788d59e960SJed Brown Input Arguments: 4798d59e960SJed Brown . adapt - time step adaptivity context 4808d59e960SJed Brown 4818d59e960SJed Brown Output Arguments: 4828d59e960SJed Brown + n - number of candidate schemes, always at least 1 4838d59e960SJed Brown . order - the order of each candidate scheme 4848d59e960SJed Brown . stageorder - the stage order of each candidate scheme 4858d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 4868d59e960SJed Brown - cost - the relative cost of each scheme 4878d59e960SJed Brown 4888d59e960SJed Brown Level: developer 4898d59e960SJed Brown 4908d59e960SJed Brown Note: 4918d59e960SJed Brown The current scheme is always returned in the first slot 4928d59e960SJed Brown 4938d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 4948d59e960SJed Brown @*/ 4958d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 4968d59e960SJed Brown { 4978d59e960SJed Brown PetscFunctionBegin; 4988d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4998d59e960SJed Brown if (n) *n = adapt->candidates.n; 5008d59e960SJed Brown if (order) *order = adapt->candidates.order; 5018d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 5028d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 5038d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 5048d59e960SJed Brown PetscFunctionReturn(0); 5058d59e960SJed Brown } 5068d59e960SJed Brown 5078d59e960SJed Brown #undef __FUNCT__ 50884df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 50984df9cb4SJed Brown /*@C 51084df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 51184df9cb4SJed Brown 51284df9cb4SJed Brown Logically Collective 51384df9cb4SJed Brown 51484df9cb4SJed Brown Input Arguments: 51584df9cb4SJed Brown + adapt - adaptive contoller 51684df9cb4SJed Brown - h - current step size 51784df9cb4SJed Brown 51884df9cb4SJed Brown Output Arguments: 51984df9cb4SJed Brown + next_sc - scheme to use for the next step 52084df9cb4SJed Brown . next_h - step size to use for the next step 52184df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 52284df9cb4SJed Brown 5231c3436cfSJed Brown Note: 5241c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 5251c3436cfSJed Brown being retried after an initial rejection. 5261c3436cfSJed Brown 52784df9cb4SJed Brown Level: developer 52884df9cb4SJed Brown 52984df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 53084df9cb4SJed Brown @*/ 53184df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 53284df9cb4SJed Brown { 53384df9cb4SJed Brown PetscErrorCode ierr; 5340b99f514SJed Brown PetscReal wlte = -1.0; 53584df9cb4SJed Brown 53684df9cb4SJed Brown PetscFunctionBegin; 53784df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 53884df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 53984df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 54084df9cb4SJed Brown PetscValidPointer(next_h,5); 54184df9cb4SJed Brown PetscValidIntPointer(accept,6); 542ce94432eSBarry Smith if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 543ce94432eSBarry 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); 5440b99f514SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 54549354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 54649354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 54749354f04SShri Abhyankar PetscReal max_time = ts->max_time; 54849354f04SShri Abhyankar PetscReal next_dt = 0.0; 54949354f04SShri Abhyankar if (ts->ptime + ts->time_step + *next_h >= max_time) { 55049354f04SShri Abhyankar next_dt = max_time - (ts->ptime + ts->time_step); 5518709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 55249354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 55349354f04SShri Abhyankar } 55449354f04SShri Abhyankar } 555ce94432eSBarry 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); 556ce94432eSBarry Smith if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h); 5571c3436cfSJed Brown 5581c3436cfSJed Brown if (adapt->monitor) { 5591c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5600b99f514SJed Brown if (wlte < 0) { 5616de24e2aSJed 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); 5620b99f514SJed Brown } else { 5636de24e2aSJed 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); 5640b99f514SJed Brown } 5651c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5661c3436cfSJed Brown } 56784df9cb4SJed Brown PetscFunctionReturn(0); 56884df9cb4SJed Brown } 56984df9cb4SJed Brown 57084df9cb4SJed Brown #undef __FUNCT__ 57197335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage" 57297335746SJed Brown /*@ 57397335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 57497335746SJed Brown 57597335746SJed Brown Collective 57697335746SJed Brown 57797335746SJed Brown Input Arguments: 57897335746SJed Brown + adapt - adaptive controller context 57997335746SJed Brown - ts - time stepper 58097335746SJed Brown 58197335746SJed Brown Output Arguments: 58297335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 58397335746SJed Brown 58497335746SJed Brown Level: developer 58597335746SJed Brown 58697335746SJed Brown .seealso: 58797335746SJed Brown @*/ 58897335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept) 58997335746SJed Brown { 59097335746SJed Brown PetscErrorCode ierr; 59197335746SJed Brown SNES snes; 59297335746SJed Brown SNESConvergedReason snesreason; 59397335746SJed Brown 59497335746SJed Brown PetscFunctionBegin; 59597335746SJed Brown *accept = PETSC_TRUE; 59697335746SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 59797335746SJed Brown ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 59897335746SJed Brown if (snesreason < 0) { 59997335746SJed Brown PetscReal dt,new_dt; 60097335746SJed Brown *accept = PETSC_FALSE; 60197335746SJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 6026de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 60397335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 60497335746SJed 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); 60597335746SJed Brown if (adapt->monitor) { 60697335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 60797335746SJed 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); 60897335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 60997335746SJed Brown } 61097335746SJed Brown } else { 61197335746SJed Brown new_dt = dt*adapt->scale_solve_failed; 61297335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 61397335746SJed Brown if (adapt->monitor) { 61497335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 61597335746SJed 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); 61697335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 61797335746SJed Brown } 61897335746SJed Brown } 61997335746SJed Brown } 6200873808bSJed Brown if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);} 62197335746SJed Brown PetscFunctionReturn(0); 62297335746SJed Brown } 62397335746SJed Brown 6240873808bSJed Brown 6250873808bSJed Brown 62697335746SJed Brown #undef __FUNCT__ 62784df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 62884df9cb4SJed Brown /*@ 62984df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 63084df9cb4SJed Brown 63184df9cb4SJed Brown Collective on MPI_Comm 63284df9cb4SJed Brown 63384df9cb4SJed Brown Input Parameter: 63484df9cb4SJed Brown . comm - The communicator 63584df9cb4SJed Brown 63684df9cb4SJed Brown Output Parameter: 63784df9cb4SJed Brown . adapt - new TSAdapt object 63884df9cb4SJed Brown 63984df9cb4SJed Brown Level: developer 64084df9cb4SJed Brown 64184df9cb4SJed Brown Notes: 64284df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 64384df9cb4SJed Brown 64484df9cb4SJed Brown .keywords: TSAdapt, create 645ad6bc421SBarry Smith .seealso: TSGetTSAdapt(), TSAdaptSetType(), TSAdaptDestroy() 64684df9cb4SJed Brown @*/ 64784df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 64884df9cb4SJed Brown { 64984df9cb4SJed Brown PetscErrorCode ierr; 65084df9cb4SJed Brown TSAdapt adapt; 65184df9cb4SJed Brown 65284df9cb4SJed Brown PetscFunctionBegin; 65384df9cb4SJed Brown *inadapt = 0; 65467c2884eSBarry Smith ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 6551c3436cfSJed Brown 6561c3436cfSJed Brown adapt->dt_min = 1e-20; 6571c3436cfSJed Brown adapt->dt_max = 1e50; 65897335746SJed Brown adapt->scale_solve_failed = 0.25; 6591c3436cfSJed Brown 66084df9cb4SJed Brown *inadapt = adapt; 66184df9cb4SJed Brown PetscFunctionReturn(0); 66284df9cb4SJed Brown } 663