184df9cb4SJed Brown 2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 384df9cb4SJed Brown 484df9cb4SJed Brown static PetscFList TSAdaptList; 584df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 684df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 784df9cb4SJed Brown static PetscClassId TSADAPT_CLASSID; 884df9cb4SJed Brown 984df9cb4SJed Brown EXTERN_C_BEGIN 1084df9cb4SJed Brown PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 11b0f836d7SJed Brown PetscErrorCode TSAdaptCreate_None(TSAdapt); 128d59e960SJed Brown PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 1384df9cb4SJed Brown EXTERN_C_END 1484df9cb4SJed Brown 1584df9cb4SJed Brown #undef __FUNCT__ 1684df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister" 1784df9cb4SJed Brown /*@C 1884df9cb4SJed Brown TSAdaptRegister - see TSAdaptRegisterDynamic() 1984df9cb4SJed Brown 2084df9cb4SJed Brown Level: advanced 2184df9cb4SJed Brown @*/ 2284df9cb4SJed Brown PetscErrorCode TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt)) 2384df9cb4SJed Brown { 2484df9cb4SJed Brown PetscErrorCode ierr; 2584df9cb4SJed Brown char fullname[PETSC_MAX_PATH_LEN]; 2684df9cb4SJed Brown 2784df9cb4SJed Brown PetscFunctionBegin; 2884df9cb4SJed Brown ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2984df9cb4SJed Brown ierr = PetscFListAdd(&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr); 3084df9cb4SJed Brown PetscFunctionReturn(0); 3184df9cb4SJed Brown } 3284df9cb4SJed Brown 3384df9cb4SJed Brown #undef __FUNCT__ 3484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll" 3584df9cb4SJed Brown /*@C 3684df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 3784df9cb4SJed Brown 3884df9cb4SJed Brown Not Collective 3984df9cb4SJed Brown 4084df9cb4SJed Brown Level: advanced 4184df9cb4SJed Brown 4284df9cb4SJed Brown .keywords: TSAdapt, register, all 4384df9cb4SJed Brown 4484df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 4584df9cb4SJed Brown @*/ 4684df9cb4SJed Brown PetscErrorCode TSAdaptRegisterAll(const char path[]) 4784df9cb4SJed Brown { 4884df9cb4SJed Brown PetscErrorCode ierr; 4984df9cb4SJed Brown 5084df9cb4SJed Brown PetscFunctionBegin; 5184df9cb4SJed Brown ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr); 52b0f836d7SJed Brown ierr = TSAdaptRegisterDynamic(TSADAPTNONE, path,"TSAdaptCreate_None", TSAdaptCreate_None);CHKERRQ(ierr); 538d59e960SJed Brown ierr = TSAdaptRegisterDynamic(TSADAPTCFL, path,"TSAdaptCreate_CFL", TSAdaptCreate_CFL);CHKERRQ(ierr); 5484df9cb4SJed Brown PetscFunctionReturn(0); 5584df9cb4SJed Brown } 5684df9cb4SJed Brown 5784df9cb4SJed Brown #undef __FUNCT__ 5884df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage" 5984df9cb4SJed Brown /*@C 6084df9cb4SJed Brown TSFinalizePackage - This function destroys everything in the TS package. It is 6184df9cb4SJed Brown called from PetscFinalize(). 6284df9cb4SJed Brown 6384df9cb4SJed Brown Level: developer 6484df9cb4SJed Brown 6584df9cb4SJed Brown .keywords: Petsc, destroy, package 6684df9cb4SJed Brown .seealso: PetscFinalize() 6784df9cb4SJed Brown @*/ 6884df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 6984df9cb4SJed Brown { 7084df9cb4SJed Brown PetscFunctionBegin; 7184df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 7284df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 7384df9cb4SJed Brown TSAdaptList = PETSC_NULL; 7484df9cb4SJed Brown PetscFunctionReturn(0); 7584df9cb4SJed Brown } 7684df9cb4SJed Brown 7784df9cb4SJed Brown #undef __FUNCT__ 7884df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage" 7984df9cb4SJed Brown /*@C 8084df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 8184df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 8284df9cb4SJed Brown TSCreate_GL() when using static libraries. 8384df9cb4SJed Brown 8484df9cb4SJed Brown Input Parameter: 8584df9cb4SJed Brown path - The dynamic library path, or PETSC_NULL 8684df9cb4SJed Brown 8784df9cb4SJed Brown Level: developer 8884df9cb4SJed Brown 8984df9cb4SJed Brown .keywords: TSAdapt, initialize, package 9084df9cb4SJed Brown .seealso: PetscInitialize() 9184df9cb4SJed Brown @*/ 9284df9cb4SJed Brown PetscErrorCode TSAdaptInitializePackage(const char path[]) 9384df9cb4SJed Brown { 9484df9cb4SJed Brown PetscErrorCode ierr; 9584df9cb4SJed Brown 9684df9cb4SJed Brown PetscFunctionBegin; 9784df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 9884df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 9984df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 10084df9cb4SJed Brown ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr); 10184df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 10284df9cb4SJed Brown PetscFunctionReturn(0); 10384df9cb4SJed Brown } 10484df9cb4SJed Brown 10584df9cb4SJed Brown #undef __FUNCT__ 10684df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterDestroy" 10784df9cb4SJed Brown /*@C 10884df9cb4SJed Brown TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic(). 10984df9cb4SJed Brown 11084df9cb4SJed Brown Not Collective 11184df9cb4SJed Brown 11284df9cb4SJed Brown Level: advanced 11384df9cb4SJed Brown 11484df9cb4SJed Brown .keywords: TSAdapt, register, destroy 11584df9cb4SJed Brown .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic() 11684df9cb4SJed Brown @*/ 11784df9cb4SJed Brown PetscErrorCode TSAdaptRegisterDestroy(void) 11884df9cb4SJed Brown { 11984df9cb4SJed Brown PetscErrorCode ierr; 12084df9cb4SJed Brown 12184df9cb4SJed Brown PetscFunctionBegin; 12284df9cb4SJed Brown ierr = PetscFListDestroy(&TSAdaptList);CHKERRQ(ierr); 12384df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 12484df9cb4SJed Brown PetscFunctionReturn(0); 12584df9cb4SJed Brown } 12684df9cb4SJed Brown 12784df9cb4SJed Brown 12884df9cb4SJed Brown #undef __FUNCT__ 12984df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 13019fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 13184df9cb4SJed Brown { 13284df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 13384df9cb4SJed Brown 13484df9cb4SJed Brown PetscFunctionBegin; 13584df9cb4SJed Brown ierr = PetscFListFind(TSAdaptList,((PetscObject)adapt)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr); 13684df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 13784df9cb4SJed Brown if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 13884df9cb4SJed Brown ierr = (*r)(adapt);CHKERRQ(ierr); 13984df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 14084df9cb4SJed Brown PetscFunctionReturn(0); 14184df9cb4SJed Brown } 14284df9cb4SJed Brown 14384df9cb4SJed Brown #undef __FUNCT__ 14484df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 14584df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 14684df9cb4SJed Brown { 14784df9cb4SJed Brown PetscErrorCode ierr; 14884df9cb4SJed Brown 14984df9cb4SJed Brown PetscFunctionBegin; 15084df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 15184df9cb4SJed Brown PetscFunctionReturn(0); 15284df9cb4SJed Brown } 15384df9cb4SJed Brown 15484df9cb4SJed Brown #undef __FUNCT__ 15584df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 15684df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 15784df9cb4SJed Brown { 15884df9cb4SJed Brown PetscErrorCode ierr; 15984df9cb4SJed Brown PetscBool iascii; 16084df9cb4SJed Brown 16184df9cb4SJed Brown PetscFunctionBegin; 162251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 16384df9cb4SJed Brown if (iascii) { 16484df9cb4SJed Brown ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr); 16584df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 16684df9cb4SJed Brown if (adapt->ops->view) { 16784df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 16884df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 16984df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 17084df9cb4SJed Brown } 17184df9cb4SJed Brown } else { 172*f2c2a1b9SBarry Smith if (adapt->ops->view) { 173*f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 174*f2c2a1b9SBarry Smith } 17584df9cb4SJed Brown } 17684df9cb4SJed Brown PetscFunctionReturn(0); 17784df9cb4SJed Brown } 17884df9cb4SJed Brown 17984df9cb4SJed Brown #undef __FUNCT__ 18084df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 18184df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 18284df9cb4SJed Brown { 18384df9cb4SJed Brown PetscErrorCode ierr; 18484df9cb4SJed Brown 18584df9cb4SJed Brown PetscFunctionBegin; 18684df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 18784df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 18884df9cb4SJed Brown if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 18984df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 1901c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 19184df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 19284df9cb4SJed Brown PetscFunctionReturn(0); 19384df9cb4SJed Brown } 19484df9cb4SJed Brown 19584df9cb4SJed Brown #undef __FUNCT__ 1961c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 1971c3436cfSJed Brown /*@ 1981c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 1991c3436cfSJed Brown 2001c3436cfSJed Brown Collective on TSAdapt 2011c3436cfSJed Brown 2021c3436cfSJed Brown Input Arguments: 2031c3436cfSJed Brown + adapt - adaptive controller context 2041c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2051c3436cfSJed Brown 2061c3436cfSJed Brown Level: intermediate 2071c3436cfSJed Brown 2081c3436cfSJed Brown .seealso: TSAdaptChoose() 2091c3436cfSJed Brown @*/ 2101c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 2111c3436cfSJed Brown { 2121c3436cfSJed Brown PetscErrorCode ierr; 2131c3436cfSJed Brown 2141c3436cfSJed Brown PetscFunctionBegin; 2151c3436cfSJed Brown if (flg) { 2161c3436cfSJed Brown if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);} 2171c3436cfSJed Brown } else { 2181c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 2191c3436cfSJed Brown } 2201c3436cfSJed Brown PetscFunctionReturn(0); 2211c3436cfSJed Brown } 2221c3436cfSJed Brown 2231c3436cfSJed Brown #undef __FUNCT__ 2240873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage" 2250873808bSJed Brown /*@C 2260873808bSJed Brown TSAdaptSetCheckStage - set a callback to check convergence for a stage 2270873808bSJed Brown 2280873808bSJed Brown Logically collective on TSAdapt 2290873808bSJed Brown 2300873808bSJed Brown Input Arguments: 2310873808bSJed Brown + adapt - adaptive controller context 2320873808bSJed Brown - func - stage check function 2330873808bSJed Brown 2340873808bSJed Brown Arguments of func: 2350873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 2360873808bSJed Brown 2370873808bSJed Brown + adapt - adaptive controller context 2380873808bSJed Brown . ts - time stepping context 2390873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 2400873808bSJed Brown 2410873808bSJed Brown Level: advanced 2420873808bSJed Brown 2430873808bSJed Brown .seealso: TSAdaptChoose() 2440873808bSJed Brown @*/ 2450873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 2460873808bSJed Brown { 2470873808bSJed Brown 2480873808bSJed Brown PetscFunctionBegin; 2490873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2500873808bSJed Brown adapt->ops->checkstage = func; 2510873808bSJed Brown PetscFunctionReturn(0); 2520873808bSJed Brown } 2530873808bSJed Brown 2540873808bSJed Brown #undef __FUNCT__ 2551c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 2561c3436cfSJed Brown /*@ 2571c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 2581c3436cfSJed Brown 2591c3436cfSJed Brown Logically Collective 2601c3436cfSJed Brown 2611c3436cfSJed Brown Input Arguments: 2621c3436cfSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 2631c3436cfSJed Brown . hmin - minimum time step 2641c3436cfSJed Brown - hmax - maximum time step 2651c3436cfSJed Brown 2661c3436cfSJed Brown Options Database Keys: 2671c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 2681c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 2691c3436cfSJed Brown 2701c3436cfSJed Brown Level: intermediate 2711c3436cfSJed Brown 2721c3436cfSJed Brown .seealso: TSAdapt 2731c3436cfSJed Brown @*/ 2741c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 2751c3436cfSJed Brown { 2761c3436cfSJed Brown 2771c3436cfSJed Brown PetscFunctionBegin; 2781c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 2791c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 2801c3436cfSJed Brown PetscFunctionReturn(0); 2811c3436cfSJed Brown } 2821c3436cfSJed Brown 2831c3436cfSJed Brown #undef __FUNCT__ 28484df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 28584df9cb4SJed Brown /*@ 28684df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 28784df9cb4SJed Brown 28884df9cb4SJed Brown Collective on TSAdapt 28984df9cb4SJed Brown 29084df9cb4SJed Brown Input Parameter: 29184df9cb4SJed Brown . adapt - the TSAdapt context 29284df9cb4SJed Brown 29384df9cb4SJed Brown Options Database Keys: 29484df9cb4SJed Brown . -ts_adapt_type <type> - basic 29584df9cb4SJed Brown 29684df9cb4SJed Brown Level: advanced 29784df9cb4SJed Brown 29884df9cb4SJed Brown Notes: 29984df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 30084df9cb4SJed Brown 30184df9cb4SJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType() 30284df9cb4SJed Brown 30384df9cb4SJed Brown .seealso: TSGetType() 30484df9cb4SJed Brown @*/ 30584df9cb4SJed Brown PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt) 30684df9cb4SJed Brown { 30784df9cb4SJed Brown PetscErrorCode ierr; 30884df9cb4SJed Brown char type[256] = TSADAPTBASIC; 3091c3436cfSJed Brown PetscBool set,flg; 31084df9cb4SJed Brown 31184df9cb4SJed Brown PetscFunctionBegin; 31284df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 31384df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 31484df9cb4SJed Brown ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr); 31584df9cb4SJed Brown ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 3168caf3d72SBarry Smith ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof(type),&flg);CHKERRQ(ierr); 31784df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 31884df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 31984df9cb4SJed Brown } 32084df9cb4SJed Brown if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);} 3211c3436cfSJed Brown ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr); 3221c3436cfSJed Brown ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr); 32397335746SJed Brown ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,PETSC_NULL);CHKERRQ(ierr); 3241c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 3251c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 32684df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 32784df9cb4SJed Brown PetscFunctionReturn(0); 32884df9cb4SJed Brown } 32984df9cb4SJed Brown 33084df9cb4SJed Brown #undef __FUNCT__ 33184df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 33284df9cb4SJed Brown /*@ 33384df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 33484df9cb4SJed Brown 33584df9cb4SJed Brown Logically Collective 33684df9cb4SJed Brown 33784df9cb4SJed Brown Input Argument: 33884df9cb4SJed Brown . adapt - adaptive controller 33984df9cb4SJed Brown 34084df9cb4SJed Brown Level: developer 34184df9cb4SJed Brown 34284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 34384df9cb4SJed Brown @*/ 34484df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 34584df9cb4SJed Brown { 34684df9cb4SJed Brown PetscErrorCode ierr; 34784df9cb4SJed Brown 34884df9cb4SJed Brown PetscFunctionBegin; 34984df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 35084df9cb4SJed Brown PetscFunctionReturn(0); 35184df9cb4SJed Brown } 35284df9cb4SJed Brown 35384df9cb4SJed Brown #undef __FUNCT__ 35484df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 35584df9cb4SJed Brown /*@C 35684df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 35784df9cb4SJed Brown 35884df9cb4SJed Brown Logically Collective 35984df9cb4SJed Brown 36084df9cb4SJed Brown Input Arguments: 36184df9cb4SJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 36284df9cb4SJed Brown . name - name of the candidate scheme to add 36384df9cb4SJed Brown . order - order of the candidate scheme 36484df9cb4SJed Brown . stageorder - stage order of the candidate scheme 3658d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 36684df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 36784df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 36884df9cb4SJed Brown 36984df9cb4SJed Brown Note: 37084df9cb4SJed Brown This routine is not available in Fortran. 37184df9cb4SJed Brown 37284df9cb4SJed Brown Level: developer 37384df9cb4SJed Brown 37484df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 37584df9cb4SJed Brown @*/ 3768d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 37784df9cb4SJed Brown { 37884df9cb4SJed Brown PetscInt c; 37984df9cb4SJed Brown 38084df9cb4SJed Brown PetscFunctionBegin; 38184df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 38284df9cb4SJed Brown if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 38384df9cb4SJed Brown if (inuse) { 38484df9cb4SJed Brown if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 38584df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 38684df9cb4SJed Brown } 3871c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 3881c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 38984df9cb4SJed Brown adapt->candidates.name[c] = name; 39084df9cb4SJed Brown adapt->candidates.order[c] = order; 39184df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 3928d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 39384df9cb4SJed Brown adapt->candidates.cost[c] = cost; 39484df9cb4SJed Brown adapt->candidates.n++; 39584df9cb4SJed Brown PetscFunctionReturn(0); 39684df9cb4SJed Brown } 39784df9cb4SJed Brown 39884df9cb4SJed Brown #undef __FUNCT__ 3998d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet" 4008d59e960SJed Brown /*@C 4018d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 4028d59e960SJed Brown 4038d59e960SJed Brown Not Collective 4048d59e960SJed Brown 4058d59e960SJed Brown Input Arguments: 4068d59e960SJed Brown . adapt - time step adaptivity context 4078d59e960SJed Brown 4088d59e960SJed Brown Output Arguments: 4098d59e960SJed Brown + n - number of candidate schemes, always at least 1 4108d59e960SJed Brown . order - the order of each candidate scheme 4118d59e960SJed Brown . stageorder - the stage order of each candidate scheme 4128d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 4138d59e960SJed Brown - cost - the relative cost of each scheme 4148d59e960SJed Brown 4158d59e960SJed Brown Level: developer 4168d59e960SJed Brown 4178d59e960SJed Brown Note: 4188d59e960SJed Brown The current scheme is always returned in the first slot 4198d59e960SJed Brown 4208d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 4218d59e960SJed Brown @*/ 4228d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 4238d59e960SJed Brown { 4248d59e960SJed Brown PetscFunctionBegin; 4258d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 4268d59e960SJed Brown if (n) *n = adapt->candidates.n; 4278d59e960SJed Brown if (order) *order = adapt->candidates.order; 4288d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 4298d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 4308d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 4318d59e960SJed Brown PetscFunctionReturn(0); 4328d59e960SJed Brown } 4338d59e960SJed Brown 4348d59e960SJed Brown #undef __FUNCT__ 43584df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 43684df9cb4SJed Brown /*@C 43784df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 43884df9cb4SJed Brown 43984df9cb4SJed Brown Logically Collective 44084df9cb4SJed Brown 44184df9cb4SJed Brown Input Arguments: 44284df9cb4SJed Brown + adapt - adaptive contoller 44384df9cb4SJed Brown - h - current step size 44484df9cb4SJed Brown 44584df9cb4SJed Brown Output Arguments: 44684df9cb4SJed Brown + next_sc - scheme to use for the next step 44784df9cb4SJed Brown . next_h - step size to use for the next step 44884df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 44984df9cb4SJed Brown 4501c3436cfSJed Brown Note: 4511c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 4521c3436cfSJed Brown being retried after an initial rejection. 4531c3436cfSJed Brown 45484df9cb4SJed Brown Level: developer 45584df9cb4SJed Brown 45684df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 45784df9cb4SJed Brown @*/ 45884df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 45984df9cb4SJed Brown { 46084df9cb4SJed Brown PetscErrorCode ierr; 4610b99f514SJed Brown PetscReal wlte = -1.0; 46284df9cb4SJed Brown 46384df9cb4SJed Brown PetscFunctionBegin; 46484df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 46584df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 46684df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 46784df9cb4SJed Brown PetscValidPointer(next_h,5); 46884df9cb4SJed Brown PetscValidIntPointer(accept,6); 46984df9cb4SJed Brown if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 47084df9cb4SJed Brown if (!adapt->candidates.inuse_set) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 4710b99f514SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 4721c3436cfSJed Brown if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1); 4731c3436cfSJed Brown if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h); 4741c3436cfSJed Brown 4751c3436cfSJed Brown if (adapt->monitor) { 4761c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 4770b99f514SJed Brown if (wlte < 0) { 4786de24e2aSJed 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); 4790b99f514SJed Brown } else { 4806de24e2aSJed 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); 4810b99f514SJed Brown } 4821c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 4831c3436cfSJed Brown } 48484df9cb4SJed Brown PetscFunctionReturn(0); 48584df9cb4SJed Brown } 48684df9cb4SJed Brown 48784df9cb4SJed Brown #undef __FUNCT__ 48897335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage" 48997335746SJed Brown /*@ 49097335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 49197335746SJed Brown 49297335746SJed Brown Collective 49397335746SJed Brown 49497335746SJed Brown Input Arguments: 49597335746SJed Brown + adapt - adaptive controller context 49697335746SJed Brown - ts - time stepper 49797335746SJed Brown 49897335746SJed Brown Output Arguments: 49997335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 50097335746SJed Brown 50197335746SJed Brown Level: developer 50297335746SJed Brown 50397335746SJed Brown .seealso: 50497335746SJed Brown @*/ 50597335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept) 50697335746SJed Brown { 50797335746SJed Brown PetscErrorCode ierr; 50897335746SJed Brown SNES snes; 50997335746SJed Brown SNESConvergedReason snesreason; 51097335746SJed Brown 51197335746SJed Brown PetscFunctionBegin; 51297335746SJed Brown *accept = PETSC_TRUE; 51397335746SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 51497335746SJed Brown ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 51597335746SJed Brown if (snesreason < 0) { 51697335746SJed Brown PetscReal dt,new_dt; 51797335746SJed Brown *accept = PETSC_FALSE; 51897335746SJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 5196de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 52097335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 52197335746SJed 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); 52297335746SJed Brown if (adapt->monitor) { 52397335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 52497335746SJed 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); 52597335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 52697335746SJed Brown } 52797335746SJed Brown } else { 52897335746SJed Brown new_dt = dt*adapt->scale_solve_failed; 52997335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 53097335746SJed Brown if (adapt->monitor) { 53197335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 53297335746SJed 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); 53397335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 53497335746SJed Brown } 53597335746SJed Brown } 53697335746SJed Brown } 5370873808bSJed Brown if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);} 53897335746SJed Brown PetscFunctionReturn(0); 53997335746SJed Brown } 54097335746SJed Brown 5410873808bSJed Brown 5420873808bSJed Brown 54397335746SJed Brown #undef __FUNCT__ 54484df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 54584df9cb4SJed Brown /*@ 54684df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 54784df9cb4SJed Brown 54884df9cb4SJed Brown Collective on MPI_Comm 54984df9cb4SJed Brown 55084df9cb4SJed Brown Input Parameter: 55184df9cb4SJed Brown . comm - The communicator 55284df9cb4SJed Brown 55384df9cb4SJed Brown Output Parameter: 55484df9cb4SJed Brown . adapt - new TSAdapt object 55584df9cb4SJed Brown 55684df9cb4SJed Brown Level: developer 55784df9cb4SJed Brown 55884df9cb4SJed Brown Notes: 55984df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 56084df9cb4SJed Brown 56184df9cb4SJed Brown .keywords: TSAdapt, create 56284df9cb4SJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 56384df9cb4SJed Brown @*/ 56484df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 56584df9cb4SJed Brown { 56684df9cb4SJed Brown PetscErrorCode ierr; 56784df9cb4SJed Brown TSAdapt adapt; 56884df9cb4SJed Brown 56984df9cb4SJed Brown PetscFunctionBegin; 57084df9cb4SJed Brown *inadapt = 0; 57184df9cb4SJed Brown ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 5721c3436cfSJed Brown 5731c3436cfSJed Brown adapt->dt_min = 1e-20; 5741c3436cfSJed Brown adapt->dt_max = 1e50; 57597335746SJed Brown adapt->scale_solve_failed = 0.25; 5761c3436cfSJed Brown 57784df9cb4SJed Brown *inadapt = adapt; 57884df9cb4SJed Brown PetscFunctionReturn(0); 57984df9cb4SJed Brown } 580