184df9cb4SJed Brown 284df9cb4SJed Brown #include <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); 1184df9cb4SJed Brown EXTERN_C_END 1284df9cb4SJed Brown 1384df9cb4SJed Brown #undef __FUNCT__ 1484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister" 1584df9cb4SJed Brown /*@C 1684df9cb4SJed Brown TSAdaptRegister - see TSAdaptRegisterDynamic() 1784df9cb4SJed Brown 1884df9cb4SJed Brown Level: advanced 1984df9cb4SJed Brown @*/ 2084df9cb4SJed Brown PetscErrorCode TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt)) 2184df9cb4SJed Brown { 2284df9cb4SJed Brown PetscErrorCode ierr; 2384df9cb4SJed Brown char fullname[PETSC_MAX_PATH_LEN]; 2484df9cb4SJed Brown 2584df9cb4SJed Brown PetscFunctionBegin; 2684df9cb4SJed Brown ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2784df9cb4SJed Brown ierr = PetscFListAdd(&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr); 2884df9cb4SJed Brown PetscFunctionReturn(0); 2984df9cb4SJed Brown } 3084df9cb4SJed Brown 3184df9cb4SJed Brown #undef __FUNCT__ 3284df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll" 3384df9cb4SJed Brown /*@C 3484df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 3584df9cb4SJed Brown 3684df9cb4SJed Brown Not Collective 3784df9cb4SJed Brown 3884df9cb4SJed Brown Level: advanced 3984df9cb4SJed Brown 4084df9cb4SJed Brown .keywords: TSAdapt, register, all 4184df9cb4SJed Brown 4284df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 4384df9cb4SJed Brown @*/ 4484df9cb4SJed Brown PetscErrorCode TSAdaptRegisterAll(const char path[]) 4584df9cb4SJed Brown { 4684df9cb4SJed Brown PetscErrorCode ierr; 4784df9cb4SJed Brown 4884df9cb4SJed Brown PetscFunctionBegin; 4984df9cb4SJed Brown ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr); 5084df9cb4SJed Brown PetscFunctionReturn(0); 5184df9cb4SJed Brown } 5284df9cb4SJed Brown 5384df9cb4SJed Brown #undef __FUNCT__ 5484df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage" 5584df9cb4SJed Brown /*@C 5684df9cb4SJed Brown TSFinalizePackage - This function destroys everything in the TS package. It is 5784df9cb4SJed Brown called from PetscFinalize(). 5884df9cb4SJed Brown 5984df9cb4SJed Brown Level: developer 6084df9cb4SJed Brown 6184df9cb4SJed Brown .keywords: Petsc, destroy, package 6284df9cb4SJed Brown .seealso: PetscFinalize() 6384df9cb4SJed Brown @*/ 6484df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 6584df9cb4SJed Brown { 6684df9cb4SJed Brown PetscFunctionBegin; 6784df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 6884df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 6984df9cb4SJed Brown TSAdaptList = PETSC_NULL; 7084df9cb4SJed Brown PetscFunctionReturn(0); 7184df9cb4SJed Brown } 7284df9cb4SJed Brown 7384df9cb4SJed Brown #undef __FUNCT__ 7484df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage" 7584df9cb4SJed Brown /*@C 7684df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 7784df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 7884df9cb4SJed Brown TSCreate_GL() when using static libraries. 7984df9cb4SJed Brown 8084df9cb4SJed Brown Input Parameter: 8184df9cb4SJed Brown path - The dynamic library path, or PETSC_NULL 8284df9cb4SJed Brown 8384df9cb4SJed Brown Level: developer 8484df9cb4SJed Brown 8584df9cb4SJed Brown .keywords: TSAdapt, initialize, package 8684df9cb4SJed Brown .seealso: PetscInitialize() 8784df9cb4SJed Brown @*/ 8884df9cb4SJed Brown PetscErrorCode TSAdaptInitializePackage(const char path[]) 8984df9cb4SJed Brown { 9084df9cb4SJed Brown PetscErrorCode ierr; 9184df9cb4SJed Brown 9284df9cb4SJed Brown PetscFunctionBegin; 9384df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 9484df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 9584df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 9684df9cb4SJed Brown ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr); 9784df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 9884df9cb4SJed Brown PetscFunctionReturn(0); 9984df9cb4SJed Brown } 10084df9cb4SJed Brown 10184df9cb4SJed Brown #undef __FUNCT__ 10284df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterDestroy" 10384df9cb4SJed Brown /*@C 10484df9cb4SJed Brown TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic(). 10584df9cb4SJed Brown 10684df9cb4SJed Brown Not Collective 10784df9cb4SJed Brown 10884df9cb4SJed Brown Level: advanced 10984df9cb4SJed Brown 11084df9cb4SJed Brown .keywords: TSAdapt, register, destroy 11184df9cb4SJed Brown .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic() 11284df9cb4SJed Brown @*/ 11384df9cb4SJed Brown PetscErrorCode TSAdaptRegisterDestroy(void) 11484df9cb4SJed Brown { 11584df9cb4SJed Brown PetscErrorCode ierr; 11684df9cb4SJed Brown 11784df9cb4SJed Brown PetscFunctionBegin; 11884df9cb4SJed Brown ierr = PetscFListDestroy(&TSAdaptList);CHKERRQ(ierr); 11984df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 12084df9cb4SJed Brown PetscFunctionReturn(0); 12184df9cb4SJed Brown } 12284df9cb4SJed Brown 12384df9cb4SJed Brown 12484df9cb4SJed Brown #undef __FUNCT__ 12584df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 12684df9cb4SJed Brown PetscErrorCode TSAdaptSetType(TSAdapt adapt,const TSAdaptType type) 12784df9cb4SJed Brown { 12884df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 12984df9cb4SJed Brown 13084df9cb4SJed Brown PetscFunctionBegin; 13184df9cb4SJed Brown ierr = PetscFListFind(TSAdaptList,((PetscObject)adapt)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr); 13284df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 13384df9cb4SJed Brown if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 13484df9cb4SJed Brown ierr = (*r)(adapt);CHKERRQ(ierr); 13584df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 13684df9cb4SJed Brown PetscFunctionReturn(0); 13784df9cb4SJed Brown } 13884df9cb4SJed Brown 13984df9cb4SJed Brown #undef __FUNCT__ 14084df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 14184df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 14284df9cb4SJed Brown { 14384df9cb4SJed Brown PetscErrorCode ierr; 14484df9cb4SJed Brown 14584df9cb4SJed Brown PetscFunctionBegin; 14684df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 14784df9cb4SJed Brown PetscFunctionReturn(0); 14884df9cb4SJed Brown } 14984df9cb4SJed Brown 15084df9cb4SJed Brown #undef __FUNCT__ 15184df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 15284df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 15384df9cb4SJed Brown { 15484df9cb4SJed Brown PetscErrorCode ierr; 15584df9cb4SJed Brown PetscBool iascii; 15684df9cb4SJed Brown 15784df9cb4SJed Brown PetscFunctionBegin; 15884df9cb4SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 15984df9cb4SJed Brown if (iascii) { 16084df9cb4SJed Brown ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr); 16184df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 16284df9cb4SJed Brown if (adapt->ops->view) { 16384df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 16484df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 16584df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 16684df9cb4SJed Brown } 16784df9cb4SJed Brown } else { 16884df9cb4SJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); 16984df9cb4SJed Brown } 17084df9cb4SJed Brown PetscFunctionReturn(0); 17184df9cb4SJed Brown } 17284df9cb4SJed Brown 17384df9cb4SJed Brown #undef __FUNCT__ 17484df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 17584df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 17684df9cb4SJed Brown { 17784df9cb4SJed Brown PetscErrorCode ierr; 17884df9cb4SJed Brown 17984df9cb4SJed Brown PetscFunctionBegin; 18084df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 18184df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 18284df9cb4SJed Brown if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 18384df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 184*1c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 18584df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 18684df9cb4SJed Brown PetscFunctionReturn(0); 18784df9cb4SJed Brown } 18884df9cb4SJed Brown 18984df9cb4SJed Brown #undef __FUNCT__ 190*1c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 191*1c3436cfSJed Brown /*@ 192*1c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 193*1c3436cfSJed Brown 194*1c3436cfSJed Brown Collective on TSAdapt 195*1c3436cfSJed Brown 196*1c3436cfSJed Brown Input Arguments: 197*1c3436cfSJed Brown + adapt - adaptive controller context 198*1c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 199*1c3436cfSJed Brown 200*1c3436cfSJed Brown Level: intermediate 201*1c3436cfSJed Brown 202*1c3436cfSJed Brown .seealso: TSAdaptChoose() 203*1c3436cfSJed Brown @*/ 204*1c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 205*1c3436cfSJed Brown { 206*1c3436cfSJed Brown PetscErrorCode ierr; 207*1c3436cfSJed Brown 208*1c3436cfSJed Brown PetscFunctionBegin; 209*1c3436cfSJed Brown if (flg) { 210*1c3436cfSJed Brown if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);} 211*1c3436cfSJed Brown } else { 212*1c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 213*1c3436cfSJed Brown } 214*1c3436cfSJed Brown PetscFunctionReturn(0); 215*1c3436cfSJed Brown } 216*1c3436cfSJed Brown 217*1c3436cfSJed Brown #undef __FUNCT__ 218*1c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 219*1c3436cfSJed Brown /*@ 220*1c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 221*1c3436cfSJed Brown 222*1c3436cfSJed Brown Logically Collective 223*1c3436cfSJed Brown 224*1c3436cfSJed Brown Input Arguments: 225*1c3436cfSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 226*1c3436cfSJed Brown . hmin - minimum time step 227*1c3436cfSJed Brown - hmax - maximum time step 228*1c3436cfSJed Brown 229*1c3436cfSJed Brown Options Database Keys: 230*1c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 231*1c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 232*1c3436cfSJed Brown 233*1c3436cfSJed Brown Level: intermediate 234*1c3436cfSJed Brown 235*1c3436cfSJed Brown .seealso: TSAdapt 236*1c3436cfSJed Brown @*/ 237*1c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 238*1c3436cfSJed Brown { 239*1c3436cfSJed Brown 240*1c3436cfSJed Brown PetscFunctionBegin; 241*1c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 242*1c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 243*1c3436cfSJed Brown PetscFunctionReturn(0); 244*1c3436cfSJed Brown } 245*1c3436cfSJed Brown 246*1c3436cfSJed Brown #undef __FUNCT__ 24784df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 24884df9cb4SJed Brown /*@ 24984df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 25084df9cb4SJed Brown 25184df9cb4SJed Brown Collective on TSAdapt 25284df9cb4SJed Brown 25384df9cb4SJed Brown Input Parameter: 25484df9cb4SJed Brown . adapt - the TSAdapt context 25584df9cb4SJed Brown 25684df9cb4SJed Brown Options Database Keys: 25784df9cb4SJed Brown . -ts_adapt_type <type> - basic 25884df9cb4SJed Brown 25984df9cb4SJed Brown Level: advanced 26084df9cb4SJed Brown 26184df9cb4SJed Brown Notes: 26284df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 26384df9cb4SJed Brown 26484df9cb4SJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType() 26584df9cb4SJed Brown 26684df9cb4SJed Brown .seealso: TSGetType() 26784df9cb4SJed Brown @*/ 26884df9cb4SJed Brown PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt) 26984df9cb4SJed Brown { 27084df9cb4SJed Brown PetscErrorCode ierr; 27184df9cb4SJed Brown char type[256] = TSADAPTBASIC; 272*1c3436cfSJed Brown PetscBool set,flg; 27384df9cb4SJed Brown 27484df9cb4SJed Brown PetscFunctionBegin; 27584df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 27684df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 27784df9cb4SJed Brown ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr); 27884df9cb4SJed Brown ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 27984df9cb4SJed Brown ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof type,&flg);CHKERRQ(ierr); 28084df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 28184df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 28284df9cb4SJed Brown } 28384df9cb4SJed Brown if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);} 284*1c3436cfSJed Brown ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr); 285*1c3436cfSJed Brown ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr); 286*1c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 287*1c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 28884df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 28984df9cb4SJed Brown PetscFunctionReturn(0); 29084df9cb4SJed Brown } 29184df9cb4SJed Brown 29284df9cb4SJed Brown #undef __FUNCT__ 29384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 29484df9cb4SJed Brown /*@ 29584df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 29684df9cb4SJed Brown 29784df9cb4SJed Brown Logically Collective 29884df9cb4SJed Brown 29984df9cb4SJed Brown Input Argument: 30084df9cb4SJed Brown . adapt - adaptive controller 30184df9cb4SJed Brown 30284df9cb4SJed Brown Level: developer 30384df9cb4SJed Brown 30484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 30584df9cb4SJed Brown @*/ 30684df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 30784df9cb4SJed Brown { 30884df9cb4SJed Brown PetscErrorCode ierr; 30984df9cb4SJed Brown 31084df9cb4SJed Brown PetscFunctionBegin; 31184df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 31284df9cb4SJed Brown PetscFunctionReturn(0); 31384df9cb4SJed Brown } 31484df9cb4SJed Brown 31584df9cb4SJed Brown #undef __FUNCT__ 31684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 31784df9cb4SJed Brown /*@C 31884df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 31984df9cb4SJed Brown 32084df9cb4SJed Brown Logically Collective 32184df9cb4SJed Brown 32284df9cb4SJed Brown Input Arguments: 32384df9cb4SJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 32484df9cb4SJed Brown . name - name of the candidate scheme to add 32584df9cb4SJed Brown . order - order of the candidate scheme 32684df9cb4SJed Brown . stageorder - stage order of the candidate scheme 32784df9cb4SJed Brown . leadingerror - leading error coefficient of the candidate scheme 32884df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 32984df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 33084df9cb4SJed Brown 33184df9cb4SJed Brown Note: 33284df9cb4SJed Brown This routine is not available in Fortran. 33384df9cb4SJed Brown 33484df9cb4SJed Brown Level: developer 33584df9cb4SJed Brown 33684df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 33784df9cb4SJed Brown @*/ 33884df9cb4SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal leadingerror,PetscReal cost,PetscBool inuse) 33984df9cb4SJed Brown { 34084df9cb4SJed Brown PetscInt c; 34184df9cb4SJed Brown 34284df9cb4SJed Brown PetscFunctionBegin; 34384df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 34484df9cb4SJed Brown if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 34584df9cb4SJed Brown if (inuse) { 34684df9cb4SJed 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()"); 34784df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 34884df9cb4SJed Brown } 349*1c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 350*1c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 35184df9cb4SJed Brown adapt->candidates.name[c] = name; 35284df9cb4SJed Brown adapt->candidates.order[c] = order; 35384df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 35484df9cb4SJed Brown adapt->candidates.leadingerror[c] = leadingerror; 35584df9cb4SJed Brown adapt->candidates.cost[c] = cost; 35684df9cb4SJed Brown adapt->candidates.n++; 35784df9cb4SJed Brown PetscFunctionReturn(0); 35884df9cb4SJed Brown } 35984df9cb4SJed Brown 36084df9cb4SJed Brown #undef __FUNCT__ 36184df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 36284df9cb4SJed Brown /*@C 36384df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 36484df9cb4SJed Brown 36584df9cb4SJed Brown Logically Collective 36684df9cb4SJed Brown 36784df9cb4SJed Brown Input Arguments: 36884df9cb4SJed Brown + adapt - adaptive contoller 36984df9cb4SJed Brown - h - current step size 37084df9cb4SJed Brown 37184df9cb4SJed Brown Output Arguments: 37284df9cb4SJed Brown + next_sc - scheme to use for the next step 37384df9cb4SJed Brown . next_h - step size to use for the next step 37484df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 37584df9cb4SJed Brown 376*1c3436cfSJed Brown Note: 377*1c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 378*1c3436cfSJed Brown being retried after an initial rejection. 379*1c3436cfSJed Brown 38084df9cb4SJed Brown Level: developer 38184df9cb4SJed Brown 38284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 38384df9cb4SJed Brown @*/ 38484df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 38584df9cb4SJed Brown { 38684df9cb4SJed Brown PetscErrorCode ierr; 38784df9cb4SJed Brown 38884df9cb4SJed Brown PetscFunctionBegin; 38984df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 39084df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 39184df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 39284df9cb4SJed Brown PetscValidPointer(next_h,5); 39384df9cb4SJed Brown PetscValidIntPointer(accept,6); 39484df9cb4SJed Brown if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 39584df9cb4SJed 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); 39684df9cb4SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept);CHKERRQ(ierr); 397*1c3436cfSJed 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); 398*1c3436cfSJed Brown if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h); 399*1c3436cfSJed Brown 400*1c3436cfSJed Brown if (adapt->monitor) { 401*1c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 402*1c3436cfSJed Brown ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %s family='%s' scheme=%D:'%s' dt=%G\n",((PetscObject)adapt)->type_name,accept?"accepted":"rejected",((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],*next_h);CHKERRQ(ierr); 403*1c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 404*1c3436cfSJed Brown } 40584df9cb4SJed Brown PetscFunctionReturn(0); 40684df9cb4SJed Brown } 40784df9cb4SJed Brown 40884df9cb4SJed Brown #undef __FUNCT__ 40984df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 41084df9cb4SJed Brown /*@ 41184df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 41284df9cb4SJed Brown 41384df9cb4SJed Brown Collective on MPI_Comm 41484df9cb4SJed Brown 41584df9cb4SJed Brown Input Parameter: 41684df9cb4SJed Brown . comm - The communicator 41784df9cb4SJed Brown 41884df9cb4SJed Brown Output Parameter: 41984df9cb4SJed Brown . adapt - new TSAdapt object 42084df9cb4SJed Brown 42184df9cb4SJed Brown Level: developer 42284df9cb4SJed Brown 42384df9cb4SJed Brown Notes: 42484df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 42584df9cb4SJed Brown 42684df9cb4SJed Brown .keywords: TSAdapt, create 42784df9cb4SJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 42884df9cb4SJed Brown @*/ 42984df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 43084df9cb4SJed Brown { 43184df9cb4SJed Brown PetscErrorCode ierr; 43284df9cb4SJed Brown TSAdapt adapt; 43384df9cb4SJed Brown 43484df9cb4SJed Brown PetscFunctionBegin; 43584df9cb4SJed Brown *inadapt = 0; 43684df9cb4SJed Brown ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 437*1c3436cfSJed Brown 438*1c3436cfSJed Brown adapt->dt_min = 1e-20; 439*1c3436cfSJed Brown adapt->dt_max = 1e50; 440*1c3436cfSJed Brown 44184df9cb4SJed Brown *inadapt = adapt; 44284df9cb4SJed Brown PetscFunctionReturn(0); 44384df9cb4SJed Brown } 444