1f26ada1bSBarry Smith /* 237f753daSBarry Smith Routines to determine options set in the options database. 3f26ada1bSBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 5ac09b921SBarry Smith 62c8e378dSBarry Smith #include <petscsys.h> 7c619b03eSJed Brown #include <petscviewertypes.h> 83a3b2205SBarry Smith 9ac09b921SBarry Smith /* SUBMANSEC = Sys */ 10ac09b921SBarry Smith 119355ec05SMatthew G. Knepley typedef enum { 129355ec05SMatthew G. Knepley PETSC_OPT_CODE, 139355ec05SMatthew G. Knepley PETSC_OPT_COMMAND_LINE, 149355ec05SMatthew G. Knepley PETSC_OPT_FILE, 159355ec05SMatthew G. Knepley PETSC_OPT_ENVIRONMENT, 169355ec05SMatthew G. Knepley NUM_PETSC_OPT_SOURCE 179355ec05SMatthew G. Knepley } PetscOptionSource; 189355ec05SMatthew G. Knepley 19c5c1f447SLisandro Dalcin #define PETSC_MAX_OPTION_NAME 512 20c5929fdfSBarry Smith typedef struct _n_PetscOptions *PetscOptions; 21c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsCreate(PetscOptions *); 22b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPush(PetscOptions); 23b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPop(void); 24c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsDestroy(PetscOptions *); 252d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsCreateDefault(void); 262d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsDestroyDefault(void); 27c5929fdfSBarry Smith 282d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsHasHelp(PetscOptions, PetscBool *); 29c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHasName(PetscOptions, const char[], const char[], PetscBool *); 30c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetBool(PetscOptions, const char[], const char[], PetscBool *, PetscBool *); 319e296098SJunchao Zhang PETSC_EXTERN PetscErrorCode PetscOptionsGetBool3(PetscOptions, const char[], const char[], PetscBool3 *, PetscBool *); 322d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetInt(PetscOptions, const char[], const char[], PetscInt *, PetscBool *); 336497c311SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetMPIInt(PetscOptions, const char[], const char[], PetscMPIInt *, PetscBool *); 342d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnum(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscBool *); 352d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEList(PetscOptions, const char[], const char[], const char *const *, PetscInt, PetscInt *, PetscBool *); 36c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetReal(PetscOptions, const char[], const char[], PetscReal *, PetscBool *); 37c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalar(PetscOptions, const char[], const char[], PetscScalar *, PetscBool *); 382d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetString(PetscOptions, const char[], const char[], char[], size_t, PetscBool *); 392d747510SLisandro Dalcin 402d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetBoolArray(PetscOptions, const char[], const char[], PetscBool[], PetscInt *, PetscBool *); 412d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnumArray(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscInt *, PetscBool *); 42c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetIntArray(PetscOptions, const char[], const char[], PetscInt[], PetscInt *, PetscBool *); 43c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetRealArray(PetscOptions, const char[], const char[], PetscReal[], PetscInt *, PetscBool *); 44c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalarArray(PetscOptions, const char[], const char[], PetscScalar[], PetscInt *, PetscBool *); 45c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetStringArray(PetscOptions, const char[], const char[], char *[], PetscInt *, PetscBool *); 463a3b2205SBarry Smith 472d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsValidKey(const char[], PetscBool *); 48c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetAlias(PetscOptions, const char[], const char[]); 49c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetValue(PetscOptions, const char[], const char[]); 50c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClearValue(PetscOptions, const char[]); 512d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPair(PetscOptions, const char[], const char[], const char *[], PetscBool *); 523a3b2205SBarry Smith 532d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetAll(PetscOptions, char *[]); 54c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsAllUsed(PetscOptions, PetscInt *); 552d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsUsed(PetscOptions, const char[], PetscBool *); 56c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeft(PetscOptions); 571ab23d95SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftGet(PetscOptions, PetscInt *, char ***, char ***); 585b191818SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftRestore(PetscOptions, PetscInt *, char ***, char ***); 59c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsView(PetscOptions, PetscViewer); 604b0e389bSBarry Smith 612d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsReject(PetscOptions, const char[], const char[], const char[]); 62c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsert(PetscOptions, int *, char ***, const char[]); 63c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertFile(MPI_Comm, PetscOptions, const char[], PetscBool); 645c23ca1cSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertFileYAML(MPI_Comm, PetscOptions, const char[], PetscBool); 65c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertString(PetscOptions, const char[]); 66080f0011SToby Isaac PETSC_EXTERN PetscErrorCode PetscOptionsInsertStringYAML(PetscOptions, const char[]); 67d06005cbSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertArgs(PetscOptions, int, char **); 68c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClear(PetscOptions); 69c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPush(PetscOptions, const char[]); 70c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPop(PetscOptions); 715d0dffe5SBarry Smith 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsGetenv(MPI_Comm, const char[], char[], size_t, PetscBool *); 732d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToBool(const char[], PetscBool *); 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToInt(const char[], PetscInt *); 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToReal(const char[], PetscReal *); 762d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToScalar(const char[], PetscScalar *); 772e8a6d31SBarry Smith 78*49abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*)(const char[], const char[], PetscOptionSource, void *), void *, PetscCtxDestroyFn *); 799355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorDefault(const char[], const char[], PetscOptionSource, void *); 80081c24baSBoyana Norris 8184761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectSetOptions(PetscObject, PetscOptions); 8284761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectGetOptions(PetscObject, PetscOptions *); 8384761bfeSJed Brown 84014dd563SJed Brown PETSC_EXTERN PetscBool PetscOptionsPublish; 85e55864a3SBarry Smith 86e55864a3SBarry Smith /* 87e55864a3SBarry Smith See manual page for PetscOptionsBegin() 884416b707SBarry Smith 894416b707SBarry Smith PetscOptionsItem and PetscOptionsItems are a single option (such as ksp_type) and a collection of such single 904416b707SBarry Smith options being handled with a PetscOptionsBegin/End() 914416b707SBarry Smith 92e55864a3SBarry Smith */ 939371c9d4SSatish Balay typedef enum { 949371c9d4SSatish Balay OPTION_INT, 959371c9d4SSatish Balay OPTION_BOOL, 969371c9d4SSatish Balay OPTION_REAL, 979371c9d4SSatish Balay OPTION_FLIST, 989371c9d4SSatish Balay OPTION_STRING, 999371c9d4SSatish Balay OPTION_REAL_ARRAY, 1009371c9d4SSatish Balay OPTION_SCALAR_ARRAY, 1019371c9d4SSatish Balay OPTION_HEAD, 1029371c9d4SSatish Balay OPTION_INT_ARRAY, 1039371c9d4SSatish Balay OPTION_ELIST, 1049371c9d4SSatish Balay OPTION_BOOL_ARRAY, 1059371c9d4SSatish Balay OPTION_STRING_ARRAY 1069371c9d4SSatish Balay } PetscOptionType; 1079355ec05SMatthew G. Knepley 1084416b707SBarry Smith typedef struct _n_PetscOptionItem *PetscOptionItem; 1094416b707SBarry Smith struct _n_PetscOptionItem { 110e55864a3SBarry Smith char *option; 111e55864a3SBarry Smith char *text; 112e55864a3SBarry Smith void *data; /* used to hold the default value and then any value it is changed to by GUI */ 113e55864a3SBarry Smith PetscFunctionList flist; /* used for available values for PetscOptionsList() */ 114e55864a3SBarry Smith const char *const *list; /* used for available values for PetscOptionsEList() */ 115e55864a3SBarry Smith char nlist; /* number of entries in list */ 116e55864a3SBarry Smith char *man; 1176497c311SBarry Smith PetscInt arraylength; /* number of entries in data in the case that it is an array (of PetscInt etc), never a giant value */ 118e55864a3SBarry Smith PetscBool set; /* the user has changed this value in the GUI */ 119e55864a3SBarry Smith PetscOptionType type; 1204416b707SBarry Smith PetscOptionItem next; 121e55864a3SBarry Smith char *pman; 122e55864a3SBarry Smith void *edata; 123e55864a3SBarry Smith }; 124e55864a3SBarry Smith 1254416b707SBarry Smith typedef struct _p_PetscOptionItems { 126e55864a3SBarry Smith PetscInt count; 1274416b707SBarry Smith PetscOptionItem next; 128e55864a3SBarry Smith char *prefix, *pprefix; 129e55864a3SBarry Smith char *title; 130e55864a3SBarry Smith MPI_Comm comm; 131e55864a3SBarry Smith PetscBool printhelp, changedmethod, alreadyprinted; 132e55864a3SBarry Smith PetscObject object; 133c5929fdfSBarry Smith PetscOptions options; 1344416b707SBarry Smith } PetscOptionItems; 13530de9b25SBarry Smith 1365f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 13794bad497SJacob Faibussowitsch extern PetscOptionItems *PetscOptionsObject; /* declare this so that the PetscOptions stubs work */ 1385f80ce2aSJacob Faibussowitsch PetscErrorCode PetscOptionsBegin(MPI_Comm, const char *, const char *, const char *); 1395f80ce2aSJacob Faibussowitsch PetscErrorCode PetscObjectOptionsBegin(PetscObject); 1405f80ce2aSJacob Faibussowitsch PetscErrorCode PetscOptionsEnd(void); 1415f80ce2aSJacob Faibussowitsch #else 14230de9b25SBarry Smith /*MC 14330de9b25SBarry Smith PetscOptionsBegin - Begins a set of queries on the options database that are related and should be 1441957e957SBarry Smith displayed on the same window of a GUI that allows the user to set the options interactively. Often one should 14516a05f60SBarry Smith use `PetscObjectOptionsBegin()` rather than this call. 14630de9b25SBarry Smith 147f2ba6396SBarry Smith Synopsis: 148aaa7dc30SBarry Smith #include <petscoptions.h> 149f2ba6396SBarry Smith PetscErrorCode PetscOptionsBegin(MPI_Comm comm,const char prefix[],const char title[],const char mansec[]) 15030de9b25SBarry Smith 151fb455bf4SPatrick Sanan Collective 15230de9b25SBarry Smith 15330de9b25SBarry Smith Input Parameters: 15430de9b25SBarry Smith + comm - communicator that shares GUI 15576280437SVaclav Hapla . prefix - options prefix for all options displayed on window (optional) 15630de9b25SBarry Smith . title - short descriptive text, for example "Krylov Solver Options" 15787497f52SBarry Smith - mansec - section of manual pages for options, for example `KSP` (optional) 15830de9b25SBarry Smith 15930de9b25SBarry Smith Level: intermediate 16030de9b25SBarry Smith 16195452b02SPatrick Sanan Notes: 162d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 163d0609cedSBarry Smith 16487497f52SBarry Smith The set of queries needs to be ended by a call to `PetscOptionsEnd()`. 165fb455bf4SPatrick Sanan 16687497f52SBarry Smith One can add subheadings with `PetscOptionsHeadBegin()`. 16730de9b25SBarry Smith 16895452b02SPatrick Sanan Developer Notes: 16916a05f60SBarry Smith `PetscOptionsPublish` is set in `PetscOptionsCheckInitial_Private()` with `-saws_options`. When `PetscOptionsPublish` is set the 17016a05f60SBarry Smith loop between `PetscOptionsBegin()` and `PetscOptionsEnd()` is run THREE times with `PetscOptionsPublishCount` of values -1,0,1. 17116a05f60SBarry Smith Otherwise the loop is run ONCE with a `PetscOptionsPublishCount` of 1. 17216a05f60SBarry Smith + \-1 - `PetscOptionsInt()` etc. just call `PetscOptionsGetInt()` etc. 17316a05f60SBarry Smith . 0 - The GUI objects are created in `PetscOptionsInt()` etc. and displayed in `PetscOptionsEnd()` and the options 17416a05f60SBarry Smith database updated with user changes; `PetscOptionsGetInt()` etc. are also called. 17516a05f60SBarry Smith - 1 - `PetscOptionsInt()` etc. again call `PetscOptionsGetInt()` etc. (possibly getting new values), in addition the help message and 176fb455bf4SPatrick Sanan default values are printed if -help was given. 17716a05f60SBarry Smith When `PetscOptionsObject.changedmethod` is set this causes `PetscOptionsPublishCount` to be reset to -2 (so in the next loop iteration it is -1) 17816a05f60SBarry Smith and the whole process is repeated. This is to handle when, for example, the `KSPType` is changed thus changing the list of 17916a05f60SBarry Smith options available so they need to be redisplayed so the user can change the. Changing `PetscOptionsObjects.changedmethod` is never 180fb455bf4SPatrick Sanan currently set. 181aee2cecaSBarry Smith 1824bb2516aSBarry Smith Fortran Note: 1834bb2516aSBarry Smith Returns ierr error code per PETSc Fortran API 1844bb2516aSBarry Smith 185db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 186db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 187db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 188db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 189c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 190db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 191db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()` 19230de9b25SBarry Smith M*/ 1939371c9d4SSatish Balay #define PetscOptionsBegin(comm, prefix, mess, sec) \ 1949371c9d4SSatish Balay do { \ 1954416b707SBarry Smith PetscOptionItems PetscOptionsObjectBase; \ 1964416b707SBarry Smith PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \ 197d0609cedSBarry Smith PetscCall(PetscMemzero(PetscOptionsObject, sizeof(*PetscOptionsObject))); \ 198e55864a3SBarry Smith for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \ 1999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBegin_Private(PetscOptionsObject, comm, prefix, mess, sec)) 20030de9b25SBarry Smith 2015fefd1ebSJed Brown /*MC 2025fefd1ebSJed Brown PetscObjectOptionsBegin - Begins a set of queries on the options database that are related and should be 2035fefd1ebSJed Brown displayed on the same window of a GUI that allows the user to set the options interactively. 2045fefd1ebSJed Brown 205f2ba6396SBarry Smith Synopsis: 206aaa7dc30SBarry Smith #include <petscoptions.h> 207f2ba6396SBarry Smith PetscErrorCode PetscObjectOptionsBegin(PetscObject obj) 2085fefd1ebSJed Brown 209c3339decSBarry Smith Collective 2105fefd1ebSJed Brown 2112fe279fdSBarry Smith Input Parameter: 2125fefd1ebSJed Brown . obj - object to set options for 2135fefd1ebSJed Brown 2145fefd1ebSJed Brown Level: intermediate 2155fefd1ebSJed Brown 21695452b02SPatrick Sanan Notes: 217d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 218d0609cedSBarry Smith 21987497f52SBarry Smith Needs to be ended by a call the `PetscOptionsEnd()` 220d0609cedSBarry Smith 22187497f52SBarry Smith Can add subheadings with `PetscOptionsHeadBegin()` 2225fefd1ebSJed Brown 223db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 224db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 225db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 226db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 227c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 228db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 229db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 2305fefd1ebSJed Brown M*/ 2319371c9d4SSatish Balay #define PetscObjectOptionsBegin(obj) \ 2329371c9d4SSatish Balay do { \ 2334416b707SBarry Smith PetscOptionItems PetscOptionsObjectBase; \ 2344416b707SBarry Smith PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \ 235c5929fdfSBarry Smith PetscOptionsObject->options = ((PetscObject)obj)->options; \ 236e55864a3SBarry Smith for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \ 237dbbe0bcdSBarry Smith PetscCall(PetscObjectOptionsBegin_Private(obj, PetscOptionsObject)) 2383194b578SJed Brown 23930de9b25SBarry Smith /*MC 24030de9b25SBarry Smith PetscOptionsEnd - Ends a set of queries on the options database that are related and should be 24130de9b25SBarry Smith displayed on the same window of a GUI that allows the user to set the options interactively. 24230de9b25SBarry Smith 243f2ba6396SBarry Smith Synopsis: 244aaa7dc30SBarry Smith #include <petscoptions.h> 245f2ba6396SBarry Smith PetscErrorCode PetscOptionsEnd(void) 24630de9b25SBarry Smith 2477cdbe19fSJose E. Roman Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()` 2487cdbe19fSJose E. Roman 24930de9b25SBarry Smith Level: intermediate 25030de9b25SBarry Smith 25195452b02SPatrick Sanan Notes: 25287497f52SBarry Smith Needs to be preceded by a call to `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` 25330de9b25SBarry Smith 254d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 255d0609cedSBarry Smith 2564bb2516aSBarry Smith Fortran Note: 2574bb2516aSBarry Smith Returns ierr error code per PETSc Fortran API 2584bb2516aSBarry Smith 259db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 260db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 261db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 262db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsHeadBegin()`, 263c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 264db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 265db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()` 26630de9b25SBarry Smith M*/ 2679371c9d4SSatish Balay #define PetscOptionsEnd() \ 2689371c9d4SSatish Balay PetscCall(PetscOptionsEnd_Private(PetscOptionsObject)); \ 2699371c9d4SSatish Balay } \ 2709371c9d4SSatish Balay } \ 2719371c9d4SSatish Balay while (0) 2725f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 27330de9b25SBarry Smith 2744416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *, MPI_Comm, const char[], const char[], const char[]); 275dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject, PetscOptionItems *); 2764416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *); 277d0609cedSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHeadBegin(PetscOptionItems *, const char[]); 27830de9b25SBarry Smith 2795f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 2809371c9d4SSatish Balay template <typename... T> 2819371c9d4SSatish Balay void PetscOptionsHeadBegin(T...); 282d0609cedSBarry Smith void PetscOptionsHeadEnd(void); 2839371c9d4SSatish Balay template <typename... T> 2849371c9d4SSatish Balay PetscErrorCode PetscOptionsEnum(T...); 2859371c9d4SSatish Balay template <typename... T> 2869371c9d4SSatish Balay PetscErrorCode PetscOptionsInt(T...); 2879371c9d4SSatish Balay template <typename... T> 2889371c9d4SSatish Balay PetscErrorCode PetscOptionsBoundedInt(T...); 2899371c9d4SSatish Balay template <typename... T> 2909371c9d4SSatish Balay PetscErrorCode PetscOptionsRangeInt(T...); 2919371c9d4SSatish Balay template <typename... T> 2929371c9d4SSatish Balay PetscErrorCode PetscOptionsReal(T...); 2939371c9d4SSatish Balay template <typename... T> 2949371c9d4SSatish Balay PetscErrorCode PetscOptionsScalar(T...); 2959371c9d4SSatish Balay template <typename... T> 2969371c9d4SSatish Balay PetscErrorCode PetscOptionsName(T...); 2979371c9d4SSatish Balay template <typename... T> 2989371c9d4SSatish Balay PetscErrorCode PetscOptionsString(T...); 2999371c9d4SSatish Balay template <typename... T> 3009371c9d4SSatish Balay PetscErrorCode PetscOptionsBool(T...); 3019371c9d4SSatish Balay template <typename... T> 3029371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(T...); 3039371c9d4SSatish Balay template <typename... T> 3049371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroup(T...); 3059371c9d4SSatish Balay template <typename... T> 3069371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(T...); 3079371c9d4SSatish Balay template <typename... T> 3089371c9d4SSatish Balay PetscErrorCode PetscOptionsFList(T...); 3099371c9d4SSatish Balay template <typename... T> 3109371c9d4SSatish Balay PetscErrorCode PetscOptionsEList(T...); 3119371c9d4SSatish Balay template <typename... T> 3129371c9d4SSatish Balay PetscErrorCode PetscOptionsRealArray(T...); 3139371c9d4SSatish Balay template <typename... T> 3149371c9d4SSatish Balay PetscErrorCode PetscOptionsScalarArray(T...); 3159371c9d4SSatish Balay template <typename... T> 3169371c9d4SSatish Balay PetscErrorCode PetscOptionsIntArray(T...); 3179371c9d4SSatish Balay template <typename... T> 3189371c9d4SSatish Balay PetscErrorCode PetscOptionsStringArray(T...); 3199371c9d4SSatish Balay template <typename... T> 3209371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolArray(T...); 3219371c9d4SSatish Balay template <typename... T> 3229371c9d4SSatish Balay PetscErrorCode PetscOptionsEnumArray(T...); 3239371c9d4SSatish Balay template <typename... T> 3249371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecated(T...); 3259371c9d4SSatish Balay template <typename... T> 3269371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecatedNoObject(T...); 3275f80ce2aSJacob Faibussowitsch #else 32830de9b25SBarry Smith /*MC 329d0609cedSBarry Smith PetscOptionsHeadBegin - Puts a heading before listing any more published options. Used, for example, 33016a05f60SBarry Smith in `KSPSetFromOptions_GMRES()`. 331d0609cedSBarry Smith 33216a05f60SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 333d0609cedSBarry Smith 334d0609cedSBarry Smith Input Parameter: 335d0609cedSBarry Smith . head - the heading text 336d0609cedSBarry Smith 33787497f52SBarry Smith Level: developer 338d0609cedSBarry Smith 339d0609cedSBarry Smith Notes: 340d0609cedSBarry Smith Handles errors directly, hence does not return an error code 341d0609cedSBarry Smith 34216a05f60SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`, and `PetscOptionsObject` created in `PetscOptionsBegin()` should be the first argument 343d0609cedSBarry Smith 34487497f52SBarry Smith Must be followed by a call to `PetscOptionsHeadEnd()` in the same function. 345d0609cedSBarry Smith 346db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 347db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 348db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 349c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 350db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 351db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 352a54bb2a9SBarry Smith M*/ 3539371c9d4SSatish Balay #define PetscOptionsHeadBegin(PetscOptionsObject, head) \ 3549371c9d4SSatish Balay do { \ 35548a46eb9SPierre Jolivet if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " %s\n", head)); \ 356d0609cedSBarry Smith } while (0) 357d0609cedSBarry Smith 358edd03b47SJacob Faibussowitsch #define PetscOptionsHead(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadBegin()", ) PetscOptionsHeadBegin(__VA_ARGS__) 359d0609cedSBarry Smith 360d0609cedSBarry Smith /*MC 36187497f52SBarry Smith PetscOptionsHeadEnd - Ends a section of options begun with `PetscOptionsHeadBegin()` 36216a05f60SBarry Smith See, for example, `KSPSetFromOptions_GMRES()`. 36330de9b25SBarry Smith 364f2ba6396SBarry Smith Synopsis: 365aaa7dc30SBarry Smith #include <petscoptions.h> 366d0609cedSBarry Smith PetscErrorCode PetscOptionsHeadEnd(void) 36730de9b25SBarry Smith 3687cdbe19fSJose E. Roman Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()` 3697cdbe19fSJose E. Roman 37030de9b25SBarry Smith Level: intermediate 37130de9b25SBarry Smith 37295452b02SPatrick Sanan Notes: 37387497f52SBarry Smith Must be between a `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` and a `PetscOptionsEnd()` 37430de9b25SBarry Smith 37587497f52SBarry Smith Must be preceded by a call to `PetscOptionsHeadBegin()` in the same function. 37630de9b25SBarry Smith 37787497f52SBarry Smith This needs to be used only if the code below `PetscOptionsHeadEnd()` can be run ONLY once. 37816a05f60SBarry Smith See, for example, `PCSetFromOptions_Composite()`. This is a `return(0)` in it for early exit 379b52f573bSBarry Smith from the function. 380b52f573bSBarry Smith 38156752e42SBarry Smith This is only for use with the PETSc options GUI 382b52f573bSBarry Smith 383db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 384db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 385db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 386c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 387db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 388db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 38930de9b25SBarry Smith M*/ 3909371c9d4SSatish Balay #define PetscOptionsHeadEnd() \ 3919371c9d4SSatish Balay do { \ 3923ba16761SJacob Faibussowitsch if (PetscOptionsObject->count != 1) PetscFunctionReturn(PETSC_SUCCESS); \ 3939371c9d4SSatish Balay } while (0) 394d0609cedSBarry Smith 395edd03b47SJacob Faibussowitsch #define PetscOptionsTail(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadEnd()", ) PetscOptionsHeadEnd(__VA_ARGS__) 396186905e3SBarry Smith 3975305534fSZach Atkins /*MC 3985305534fSZach Atkins PetscOptionsEnum - Gets the enum value for a particular option in the database. 3995305534fSZach Atkins 4005305534fSZach Atkins Synopsis: 4015305534fSZach Atkins #include <petscoptions.h> 4025305534fSZach Atkins PetscErrorCode PetscOptionsEnum(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum currentvalue, PetscEnum *value, PetscBool *set) 4035305534fSZach Atkins 4045305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 4055305534fSZach Atkins 4065305534fSZach Atkins Input Parameters: 4075305534fSZach Atkins + opt - option name 4085305534fSZach Atkins . text - short string that describes the option 4095305534fSZach Atkins . man - manual page with additional information on option 4105305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 4115305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 4125305534fSZach Atkins .vb 4135305534fSZach Atkins PetscOptionsEnum(..., obj->value,&object->value,...) or 4145305534fSZach Atkins value = defaultvalue 4159314d9b7SBarry Smith PetscOptionsEnum(..., value,&value,&set); 4169314d9b7SBarry Smith if (set) { 4175305534fSZach Atkins .ve 4185305534fSZach Atkins 4195305534fSZach Atkins Output Parameters: 4205305534fSZach Atkins + value - the value to return 4215305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 4225305534fSZach Atkins 4235305534fSZach Atkins Level: beginner 4245305534fSZach Atkins 4255305534fSZach Atkins Notes: 4265305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 4275305534fSZach Atkins 4289314d9b7SBarry Smith `list` is usually something like `PCASMTypes` or some other predefined list of enum names 4295305534fSZach Atkins 4305305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 4319314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 4325305534fSZach Atkins 4335305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 4345305534fSZach Atkins 4355305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 4365305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 4375305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 4385305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 4395305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 4405305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 4415305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 4425305534fSZach Atkins M*/ 4431ff8fb82SZach Atkins #define PetscOptionsEnum(opt, text, man, list, currentvalue, value, set) PetscOptionsEnum_Private(PetscOptionsObject, opt, text, man, list, currentvalue, value, set) 4445305534fSZach Atkins 4455305534fSZach Atkins /*MC 4465305534fSZach Atkins PetscOptionsInt - Gets the integer value for a particular option in the database. 4475305534fSZach Atkins 4485305534fSZach Atkins Synopsis: 4495305534fSZach Atkins #include <petscoptions.h> 4505305534fSZach Atkins PetscErrorCode PetscOptionsInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set) 4515305534fSZach Atkins 4525305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 4535305534fSZach Atkins 4545305534fSZach Atkins Input Parameters: 4555305534fSZach Atkins + opt - option name 4565305534fSZach Atkins . text - short string that describes the option 4575305534fSZach Atkins . man - manual page with additional information on option 4585305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 4595305534fSZach Atkins .vb 4605305534fSZach Atkins PetscOptionsInt(..., obj->value, &obj->value, ...) or 4615305534fSZach Atkins value = defaultvalue 4629314d9b7SBarry Smith PetscOptionsInt(..., value, &value, &set); 4639314d9b7SBarry Smith if (set) { 4645305534fSZach Atkins .ve 4655305534fSZach Atkins 4665305534fSZach Atkins Output Parameters: 4675305534fSZach Atkins + value - the integer value to return 4685305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 4695305534fSZach Atkins 4705305534fSZach Atkins Level: beginner 4715305534fSZach Atkins 4725305534fSZach Atkins Notes: 4735305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 4749314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 4755305534fSZach Atkins 4765305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 4775305534fSZach Atkins 4785305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 4795305534fSZach Atkins 4805305534fSZach Atkins .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 4815305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 4825305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 4835305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 4845305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 4855305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 48652ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 4875305534fSZach Atkins M*/ 4881690c2aeSBarry Smith #define PetscOptionsInt(opt, text, man, currentvalue, value, set) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_INT_MIN, PETSC_INT_MAX) 4895305534fSZach Atkins 4905305534fSZach Atkins /*MC 4916497c311SBarry Smith PetscOptionsMPIInt - Gets the MPI integer value for a particular option in the database. 4926497c311SBarry Smith 4936497c311SBarry Smith Synopsis: 4946497c311SBarry Smith #include <petscoptions.h> 4956497c311SBarry Smith PetscErrorCode PetscOptionsMPIInt(const char opt[], const char text[], const char man[], PetscMPIInt currentvalue, PetscMPIInt *value, PetscBool *set) 4966497c311SBarry Smith 4976497c311SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 4986497c311SBarry Smith 4996497c311SBarry Smith Input Parameters: 5006497c311SBarry Smith + opt - option name 5016497c311SBarry Smith . text - short string that describes the option 5026497c311SBarry Smith . man - manual page with additional information on option 5036497c311SBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 5046497c311SBarry Smith .vb 5056497c311SBarry Smith PetscOptionsInt(..., obj->value, &obj->value, ...) or 5066497c311SBarry Smith value = defaultvalue 5076497c311SBarry Smith PetscOptionsInt(..., value, &value, &set); 5086497c311SBarry Smith if (set) { 5096497c311SBarry Smith .ve 5106497c311SBarry Smith 5116497c311SBarry Smith Output Parameters: 5126497c311SBarry Smith + value - the MPI integer value to return 5136497c311SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 5146497c311SBarry Smith 5156497c311SBarry Smith Level: beginner 5166497c311SBarry Smith 5176497c311SBarry Smith Notes: 5186497c311SBarry Smith If the user does not supply the option at all `value` is NOT changed. Thus 5196497c311SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 5206497c311SBarry Smith 5216497c311SBarry Smith The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 5226497c311SBarry Smith 5236497c311SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 5246497c311SBarry Smith 5256497c311SBarry Smith .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 5266497c311SBarry Smith `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 5276497c311SBarry Smith `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 5286497c311SBarry Smith `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 5296497c311SBarry Smith `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 5306497c311SBarry Smith `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 5316497c311SBarry Smith `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 5326497c311SBarry Smith M*/ 5336497c311SBarry Smith #define PetscOptionsMPIInt(opt, text, man, currentvalue, value, set) PetscOptionsMPIInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MPI_INT_MIN, PETSC_MPI_INT_MAX) 5346497c311SBarry Smith 5356497c311SBarry Smith /*MC 53652ce0ab5SPierre Jolivet PetscOptionsBoundedInt - Gets an integer value greater than or equal to a given bound for a particular option in the database. 5375305534fSZach Atkins 5385305534fSZach Atkins Synopsis: 5395305534fSZach Atkins #include <petscoptions.h> 5409314d9b7SBarry Smith PetscErrorCode PetscOptionsBoundedInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt bound) 5415305534fSZach Atkins 5425305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 5435305534fSZach Atkins 5445305534fSZach Atkins Input Parameters: 5455305534fSZach Atkins + opt - option name 5465305534fSZach Atkins . text - short string that describes the option 5475305534fSZach Atkins . man - manual page with additional information on option 5485305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 5495305534fSZach Atkins .vb 55052ce0ab5SPierre Jolivet PetscOptionsBoundedInt(..., obj->value, &obj->value, ...) 5515305534fSZach Atkins .ve 5525305534fSZach Atkins or 5535305534fSZach Atkins .vb 5545305534fSZach Atkins value = defaultvalue 55552ce0ab5SPierre Jolivet PetscOptionsBoundedInt(..., value, &value, &set, ...); 5569314d9b7SBarry Smith if (set) { 5575305534fSZach Atkins .ve 55852ce0ab5SPierre Jolivet - bound - the requested value should be greater than or equal to this bound or an error is generated 5595305534fSZach Atkins 5605305534fSZach Atkins Output Parameters: 5615305534fSZach Atkins + value - the integer value to return 5629314d9b7SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 5635305534fSZach Atkins 5645305534fSZach Atkins Level: beginner 5655305534fSZach Atkins 5665305534fSZach Atkins Notes: 5675305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 5689314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 5695305534fSZach Atkins 5705305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 5715305534fSZach Atkins 572af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 5735305534fSZach Atkins 5745305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 5755305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 5765305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 5775305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 5785305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 5795305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 58052ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 5815305534fSZach Atkins M*/ 5821690c2aeSBarry Smith #define PetscOptionsBoundedInt(opt, text, man, currentvalue, value, set, lb) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_INT_MAX) 5835305534fSZach Atkins 5845305534fSZach Atkins /*MC 5855305534fSZach Atkins PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database. 5865305534fSZach Atkins 5875305534fSZach Atkins Synopsis: 5885305534fSZach Atkins #include <petscoptions.h> 5899314d9b7SBarry Smith PetscErrorCode PetscOptionsRangeInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt lb, PetscInt ub) 5905305534fSZach Atkins 5915305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 5925305534fSZach Atkins 5935305534fSZach Atkins Input Parameters: 5945305534fSZach Atkins + opt - option name 5955305534fSZach Atkins . text - short string that describes the option 5965305534fSZach Atkins . man - manual page with additional information on option 5975305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 5985305534fSZach Atkins .vb 59952ce0ab5SPierre Jolivet PetscOptionsRangeInt(..., obj->value, &obj->value, ...) 60052ce0ab5SPierre Jolivet .ve 60152ce0ab5SPierre Jolivet or 60252ce0ab5SPierre Jolivet .vb 6035305534fSZach Atkins value = defaultvalue 60452ce0ab5SPierre Jolivet PetscOptionsRangeInt(..., value, &value, &set, ...); 6059314d9b7SBarry Smith if (set) { 6065305534fSZach Atkins .ve 6075305534fSZach Atkins . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 6085305534fSZach Atkins - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 6095305534fSZach Atkins 6105305534fSZach Atkins Output Parameters: 6115305534fSZach Atkins + value - the integer value to return 6129314d9b7SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 6135305534fSZach Atkins 6145305534fSZach Atkins Level: beginner 6155305534fSZach Atkins 6165305534fSZach Atkins Notes: 6175305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 6189314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 6195305534fSZach Atkins 6205305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 6215305534fSZach Atkins 622af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 6235305534fSZach Atkins 6245305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 6255305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 6265305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 6275305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 6285305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 6295305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 63052ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 6315305534fSZach Atkins M*/ 6321ff8fb82SZach Atkins #define PetscOptionsRangeInt(opt, text, man, currentvalue, value, set, lb, ub) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub) 6335305534fSZach Atkins 6345305534fSZach Atkins /*MC 63552ce0ab5SPierre Jolivet PetscOptionsReal - Gets a `PetscReal` value for a particular option in the database. 6365305534fSZach Atkins 6375305534fSZach Atkins Synopsis: 6385305534fSZach Atkins #include <petscoptions.h> 6395305534fSZach Atkins PetscErrorCode PetscOptionsReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set) 6405305534fSZach Atkins 6415305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 6425305534fSZach Atkins 6435305534fSZach Atkins Input Parameters: 6445305534fSZach Atkins + opt - option name 6455305534fSZach Atkins . text - short string that describes the option 6465305534fSZach Atkins . man - manual page with additional information on option 6475305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6485305534fSZach Atkins .vb 6495305534fSZach Atkins PetscOptionsReal(..., obj->value,&obj->value,...) or 6505305534fSZach Atkins value = defaultvalue 6519314d9b7SBarry Smith PetscOptionsReal(..., value,&value,&set); 6529314d9b7SBarry Smith if (set) { 6535305534fSZach Atkins .ve 6545305534fSZach Atkins 6555305534fSZach Atkins Output Parameters: 6565305534fSZach Atkins + value - the value to return 6575305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 6585305534fSZach Atkins 6595305534fSZach Atkins Level: beginner 6605305534fSZach Atkins 6615305534fSZach Atkins Notes: 6625305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 6639314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 6645305534fSZach Atkins 6655305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 6665305534fSZach Atkins 667af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 6685305534fSZach Atkins 6695305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 6705305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 6715305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 6725305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 6735305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 6745305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 67552ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 6765305534fSZach Atkins M*/ 67752ce0ab5SPierre Jolivet #define PetscOptionsReal(opt, text, man, currentvalue, value, set) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MIN_REAL, PETSC_MAX_REAL) 67852ce0ab5SPierre Jolivet 67952ce0ab5SPierre Jolivet /*MC 68052ce0ab5SPierre Jolivet PetscOptionsBoundedReal - Gets a `PetscReal` value greater than or equal to a given bound for a particular option in the database. 68152ce0ab5SPierre Jolivet 68252ce0ab5SPierre Jolivet Synopsis: 68352ce0ab5SPierre Jolivet #include <petscoptions.h> 68452ce0ab5SPierre Jolivet PetscErrorCode PetscOptionsBoundedReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal bound) 68552ce0ab5SPierre Jolivet 68652ce0ab5SPierre Jolivet Logically Collective on the communicator passed in `PetscOptionsBegin()` 68752ce0ab5SPierre Jolivet 68852ce0ab5SPierre Jolivet Input Parameters: 68952ce0ab5SPierre Jolivet + opt - option name 69052ce0ab5SPierre Jolivet . text - short string that describes the option 69152ce0ab5SPierre Jolivet . man - manual page with additional information on option 69252ce0ab5SPierre Jolivet . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 69352ce0ab5SPierre Jolivet .vb 69452ce0ab5SPierre Jolivet PetscOptionsBoundedReal(..., obj->value, &obj->value, ...) 69552ce0ab5SPierre Jolivet .ve 69652ce0ab5SPierre Jolivet or 69752ce0ab5SPierre Jolivet .vb 69852ce0ab5SPierre Jolivet value = defaultvalue 69952ce0ab5SPierre Jolivet PetscOptionsBoundedReal(..., value, &value, &set, ...); 70052ce0ab5SPierre Jolivet if (set) { 70152ce0ab5SPierre Jolivet .ve 70252ce0ab5SPierre Jolivet - bound - the requested value should be greater than or equal to this bound or an error is generated 70352ce0ab5SPierre Jolivet 70452ce0ab5SPierre Jolivet Output Parameters: 70552ce0ab5SPierre Jolivet + value - the real value to return 70652ce0ab5SPierre Jolivet - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 70752ce0ab5SPierre Jolivet 70852ce0ab5SPierre Jolivet Level: beginner 70952ce0ab5SPierre Jolivet 71052ce0ab5SPierre Jolivet Notes: 71152ce0ab5SPierre Jolivet If the user does not supply the option at all `value` is NOT changed. Thus 71252ce0ab5SPierre Jolivet you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 71352ce0ab5SPierre Jolivet 71452ce0ab5SPierre Jolivet The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 71552ce0ab5SPierre Jolivet 71652ce0ab5SPierre Jolivet Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 71752ce0ab5SPierre Jolivet 71852ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 71952ce0ab5SPierre Jolivet `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 72052ce0ab5SPierre Jolivet `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 72152ce0ab5SPierre Jolivet `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 72252ce0ab5SPierre Jolivet `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 72352ce0ab5SPierre Jolivet `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 72452ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedInt()`, `PetscOptionsRangeReal()` 72552ce0ab5SPierre Jolivet M*/ 72652ce0ab5SPierre Jolivet #define PetscOptionsBoundedReal(opt, text, man, currentvalue, value, set, lb) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_MAX_REAL) 72752ce0ab5SPierre Jolivet 72852ce0ab5SPierre Jolivet /*MC 72952ce0ab5SPierre Jolivet PetscOptionsRangeReal - Gets a `PetscReal` value within a range of values for a particular option in the database. 73052ce0ab5SPierre Jolivet 73152ce0ab5SPierre Jolivet Synopsis: 73252ce0ab5SPierre Jolivet #include <petscoptions.h> 73352ce0ab5SPierre Jolivet PetscErrorCode PetscOptionsRangeReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal lb, PetscReal ub) 73452ce0ab5SPierre Jolivet 73552ce0ab5SPierre Jolivet Logically Collective on the communicator passed in `PetscOptionsBegin()` 73652ce0ab5SPierre Jolivet 73752ce0ab5SPierre Jolivet Input Parameters: 73852ce0ab5SPierre Jolivet + opt - option name 73952ce0ab5SPierre Jolivet . text - short string that describes the option 74052ce0ab5SPierre Jolivet . man - manual page with additional information on option 74152ce0ab5SPierre Jolivet . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 74252ce0ab5SPierre Jolivet .vb 74352ce0ab5SPierre Jolivet PetscOptionsRangeReal(..., obj->value, &obj->value, ...) 74452ce0ab5SPierre Jolivet .ve 74552ce0ab5SPierre Jolivet or 74652ce0ab5SPierre Jolivet .vb 74752ce0ab5SPierre Jolivet value = defaultvalue 74852ce0ab5SPierre Jolivet PetscOptionsRangeReal(..., value, &value, &set, ...); 74952ce0ab5SPierre Jolivet if (set) { 75052ce0ab5SPierre Jolivet .ve 75152ce0ab5SPierre Jolivet . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 75252ce0ab5SPierre Jolivet - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 75352ce0ab5SPierre Jolivet 75452ce0ab5SPierre Jolivet Output Parameters: 75552ce0ab5SPierre Jolivet + value - the value to return 75652ce0ab5SPierre Jolivet - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 75752ce0ab5SPierre Jolivet 75852ce0ab5SPierre Jolivet Level: beginner 75952ce0ab5SPierre Jolivet 76052ce0ab5SPierre Jolivet Notes: 76152ce0ab5SPierre Jolivet If the user does not supply the option at all `value` is NOT changed. Thus 76252ce0ab5SPierre Jolivet you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 76352ce0ab5SPierre Jolivet 76452ce0ab5SPierre Jolivet The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 76552ce0ab5SPierre Jolivet 76652ce0ab5SPierre Jolivet Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 76752ce0ab5SPierre Jolivet 76852ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 76952ce0ab5SPierre Jolivet `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 77052ce0ab5SPierre Jolivet `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 77152ce0ab5SPierre Jolivet `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 77252ce0ab5SPierre Jolivet `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 77352ce0ab5SPierre Jolivet `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 77452ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRangeInt()`, `PetscOptionsBoundedReal()` 77552ce0ab5SPierre Jolivet M*/ 77652ce0ab5SPierre Jolivet #define PetscOptionsRangeReal(opt, text, man, currentvalue, value, set, lb, ub) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub) 7775305534fSZach Atkins 7785305534fSZach Atkins /*MC 7795305534fSZach Atkins PetscOptionsScalar - Gets the `PetscScalar` value for a particular option in the database. 7805305534fSZach Atkins 7815305534fSZach Atkins Synopsis: 7825305534fSZach Atkins #include <petscoptions.h> 7835305534fSZach Atkins PetscErrorCode PetscOptionsScalar(const char opt[], const char text[], const char man[], PetscScalar currentvalue, PetscScalar *value, PetscBool *set) 7845305534fSZach Atkins 7855305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 7865305534fSZach Atkins 7875305534fSZach Atkins Input Parameters: 7885305534fSZach Atkins + opt - option name 7895305534fSZach Atkins . text - short string that describes the option 7905305534fSZach Atkins . man - manual page with additional information on option 7915305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 7925305534fSZach Atkins .vb 7935305534fSZach Atkins PetscOptionsScalar(..., obj->value,&obj->value,...) or 7945305534fSZach Atkins value = defaultvalue 7959314d9b7SBarry Smith PetscOptionsScalar(..., value,&value,&set); 7969314d9b7SBarry Smith if (set) { 7975305534fSZach Atkins .ve 7985305534fSZach Atkins 7995305534fSZach Atkins Output Parameters: 8005305534fSZach Atkins + value - the value to return 8015305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 8025305534fSZach Atkins 8035305534fSZach Atkins Level: beginner 8045305534fSZach Atkins 8055305534fSZach Atkins Notes: 8065305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 8079314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 8085305534fSZach Atkins 8095305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 8105305534fSZach Atkins 811af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8125305534fSZach Atkins 8135305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8145305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8155305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 8165305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8175305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8185305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8195305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8205305534fSZach Atkins M*/ 8211ff8fb82SZach Atkins #define PetscOptionsScalar(opt, text, man, currentvalue, value, set) PetscOptionsScalar_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 8225305534fSZach Atkins 8235305534fSZach Atkins /*MC 8245305534fSZach Atkins PetscOptionsName - Determines if a particular option has been set in the database. This returns true whether the option is a number, string or boolean, even 8255305534fSZach Atkins its value is set to false. 8265305534fSZach Atkins 8275305534fSZach Atkins Synopsis: 8285305534fSZach Atkins #include <petscoptions.h> 8299314d9b7SBarry Smith PetscErrorCode PetscOptionsName(const char opt[], const char text[], const char man[], PetscBool *set) 8305305534fSZach Atkins 8315305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 8325305534fSZach Atkins 8335305534fSZach Atkins Input Parameters: 8345305534fSZach Atkins + opt - option name 8355305534fSZach Atkins . text - short string that describes the option 8365305534fSZach Atkins - man - manual page with additional information on option 8375305534fSZach Atkins 8385305534fSZach Atkins Output Parameter: 8399314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 8405305534fSZach Atkins 8415305534fSZach Atkins Level: beginner 8425305534fSZach Atkins 8435305534fSZach Atkins Note: 844af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8455305534fSZach Atkins 8465305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8475305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8485305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 8495305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8505305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8515305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8525305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8535305534fSZach Atkins M*/ 8549314d9b7SBarry Smith #define PetscOptionsName(opt, text, man, set) PetscOptionsName_Private(PetscOptionsObject, opt, text, man, set) 8555305534fSZach Atkins 8565305534fSZach Atkins /*MC 8575305534fSZach Atkins PetscOptionsString - Gets the string value for a particular option in the database. 8585305534fSZach Atkins 8595305534fSZach Atkins Synopsis: 8605305534fSZach Atkins #include <petscoptions.h> 8615305534fSZach Atkins PetscErrorCode PetscOptionsString(const char opt[], const char text[], const char man[], const char currentvalue[], char value[], size_t len, PetscBool *set) 8625305534fSZach Atkins 8635305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 8645305534fSZach Atkins 8655305534fSZach Atkins Input Parameters: 8665305534fSZach Atkins + opt - option name 8675305534fSZach Atkins . text - short string that describes the option 8685305534fSZach Atkins . man - manual page with additional information on option 8695305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 8705305534fSZach Atkins - len - length of the result string including null terminator 8715305534fSZach Atkins 8725305534fSZach Atkins Output Parameters: 8735305534fSZach Atkins + value - the value to return 8745305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 8755305534fSZach Atkins 8765305534fSZach Atkins Level: beginner 8775305534fSZach Atkins 8785305534fSZach Atkins Notes: 879af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8805305534fSZach Atkins 8819314d9b7SBarry Smith If the user provided no string (for example `-optionname` `-someotheroption`) `set` is set to `PETSC_TRUE` (and the string is filled with nulls). 8825305534fSZach Atkins 8835305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 8849314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 8855305534fSZach Atkins 8865305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 8875305534fSZach Atkins 8885305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8895305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8905305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 8915305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8925305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8935305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8945305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8955305534fSZach Atkins M*/ 8961ff8fb82SZach Atkins #define PetscOptionsString(opt, text, man, currentvalue, value, len, set) PetscOptionsString_Private(PetscOptionsObject, opt, text, man, currentvalue, value, len, set) 8975305534fSZach Atkins 8985305534fSZach Atkins /*MC 8995305534fSZach Atkins PetscOptionsBool - Determines if a particular option is in the database with a true or false 9005305534fSZach Atkins 9015305534fSZach Atkins Synopsis: 9025305534fSZach Atkins #include <petscoptions.h> 9035305534fSZach Atkins PetscErrorCode PetscOptionsBool(const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool *flg, PetscBool *set) 9045305534fSZach Atkins 9055305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9065305534fSZach Atkins 9075305534fSZach Atkins Input Parameters: 9085305534fSZach Atkins + opt - option name 9095305534fSZach Atkins . text - short string that describes the option 9105305534fSZach Atkins . man - manual page with additional information on option 9115305534fSZach Atkins - currentvalue - the current value 9125305534fSZach Atkins 9135305534fSZach Atkins Output Parameters: 9145305534fSZach Atkins + flg - `PETSC_TRUE` or `PETSC_FALSE` 9155305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 9165305534fSZach Atkins 9175305534fSZach Atkins Level: beginner 9185305534fSZach Atkins 9195305534fSZach Atkins Notes: 9205305534fSZach Atkins TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE` 9215305534fSZach Atkins FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE` 9225305534fSZach Atkins 9239314d9b7SBarry Smith If the option is given, but no value is provided, then `flg` and `set` are both given the value `PETSC_TRUE`. That is `-requested_bool` 9245305534fSZach Atkins is equivalent to `-requested_bool true` 9255305534fSZach Atkins 9265305534fSZach Atkins If the user does not supply the option at all `flg` is NOT changed. Thus 9279314d9b7SBarry Smith you should ALWAYS initialize the `flg` variable if you access it without first checking that the `set` flag is `PETSC_TRUE`. 9285305534fSZach Atkins 9295305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9305305534fSZach Atkins 9315305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 9325305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 9335305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 9345305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 9355305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 9365305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 9375305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 9385305534fSZach Atkins M*/ 9391ff8fb82SZach Atkins #define PetscOptionsBool(opt, text, man, currentvalue, value, set) PetscOptionsBool_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 9405305534fSZach Atkins 9415305534fSZach Atkins /*MC 9425305534fSZach Atkins PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 9435305534fSZach Atkins which at most a single value can be true. 9445305534fSZach Atkins 9455305534fSZach Atkins Synopsis: 9465305534fSZach Atkins #include <petscoptions.h> 9479314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[], const char text[], const char man[], PetscBool *set) 9485305534fSZach Atkins 9495305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9505305534fSZach Atkins 9515305534fSZach Atkins Input Parameters: 9525305534fSZach Atkins + opt - option name 9535305534fSZach Atkins . text - short string that describes the option 9545305534fSZach Atkins - man - manual page with additional information on option 9555305534fSZach Atkins 9565305534fSZach Atkins Output Parameter: 9579314d9b7SBarry Smith . set - whether that option was set or not 9585305534fSZach Atkins 9595305534fSZach Atkins Level: intermediate 9605305534fSZach Atkins 9615305534fSZach Atkins Notes: 962af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9635305534fSZach Atkins 9645305534fSZach Atkins Must be followed by 0 or more `PetscOptionsBoolGroup()`s and `PetscOptionsBoolGroupEnd()` 9655305534fSZach Atkins 9665305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 9675305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 9685305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 9695305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 9705305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 9715305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 9725305534fSZach Atkins M*/ 9739314d9b7SBarry Smith #define PetscOptionsBoolGroupBegin(opt, text, man, set) PetscOptionsBoolGroupBegin_Private(PetscOptionsObject, opt, text, man, set) 9745305534fSZach Atkins 9755305534fSZach Atkins /*MC 9765305534fSZach Atkins PetscOptionsBoolGroup - One in a series of logical queries on the options database for 9775305534fSZach Atkins which at most a single value can be true. 9785305534fSZach Atkins 9795305534fSZach Atkins Synopsis: 9805305534fSZach Atkins #include <petscoptions.h> 9819314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroup(const char opt[], const char text[], const char man[], PetscBool *set) 9825305534fSZach Atkins 9835305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9845305534fSZach Atkins 9855305534fSZach Atkins Input Parameters: 9865305534fSZach Atkins + opt - option name 9875305534fSZach Atkins . text - short string that describes the option 9885305534fSZach Atkins - man - manual page with additional information on option 9895305534fSZach Atkins 9905305534fSZach Atkins Output Parameter: 9919314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 9925305534fSZach Atkins 9935305534fSZach Atkins Level: intermediate 9945305534fSZach Atkins 9955305534fSZach Atkins Notes: 996af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9975305534fSZach Atkins 9985305534fSZach Atkins Must follow a `PetscOptionsBoolGroupBegin()` and preceded a `PetscOptionsBoolGroupEnd()` 9995305534fSZach Atkins 10005305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 10015305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 10025305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 10035305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 10045305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 10055305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 10065305534fSZach Atkins M*/ 10079314d9b7SBarry Smith #define PetscOptionsBoolGroup(opt, text, man, set) PetscOptionsBoolGroup_Private(PetscOptionsObject, opt, text, man, set) 10085305534fSZach Atkins 10095305534fSZach Atkins /*MC 10105305534fSZach Atkins PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 10115305534fSZach Atkins which at most a single value can be true. 10125305534fSZach Atkins 10135305534fSZach Atkins Synopsis: 10145305534fSZach Atkins #include <petscoptions.h> 10159314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[], const char text[], const char man[], PetscBool *set) 10165305534fSZach Atkins 10175305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 10185305534fSZach Atkins 10195305534fSZach Atkins Input Parameters: 10205305534fSZach Atkins + opt - option name 10215305534fSZach Atkins . text - short string that describes the option 10225305534fSZach Atkins - man - manual page with additional information on option 10235305534fSZach Atkins 10245305534fSZach Atkins Output Parameter: 10259314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 10265305534fSZach Atkins 10275305534fSZach Atkins Level: intermediate 10285305534fSZach Atkins 10295305534fSZach Atkins Notes: 1030af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 10315305534fSZach Atkins 10325305534fSZach Atkins Must follow a `PetscOptionsBoolGroupBegin()` 10335305534fSZach Atkins 10345305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 10355305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 10365305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 10375305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 10385305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 10395305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 10405305534fSZach Atkins M*/ 10419314d9b7SBarry Smith #define PetscOptionsBoolGroupEnd(opt, text, man, set) PetscOptionsBoolGroupEnd_Private(PetscOptionsObject, opt, text, man, set) 10425305534fSZach Atkins 10435305534fSZach Atkins /*MC 10445305534fSZach Atkins PetscOptionsFList - Puts a list of option values that a single one may be selected from 10455305534fSZach Atkins 10465305534fSZach Atkins Synopsis: 10475305534fSZach Atkins #include <petscoptions.h> 10485305534fSZach Atkins PetscErrorCode PetscOptionsFList(const char opt[], const char ltext[], const char man[], PetscFunctionList list, const char currentvalue[], char value[], size_t len, PetscBool *set) 10495305534fSZach Atkins 10505305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 10515305534fSZach Atkins 10525305534fSZach Atkins Input Parameters: 10535305534fSZach Atkins + opt - option name 10545305534fSZach Atkins . ltext - short string that describes the option 10555305534fSZach Atkins . man - manual page with additional information on option 10565305534fSZach Atkins . list - the possible choices 10575305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10585305534fSZach Atkins .vb 10599314d9b7SBarry Smith PetscOptionsFlist(..., obj->value,value,len,&set); 10609314d9b7SBarry Smith if (set) { 10615305534fSZach Atkins .ve 10625305534fSZach Atkins - len - the length of the character array value 10635305534fSZach Atkins 10645305534fSZach Atkins Output Parameters: 10655305534fSZach Atkins + value - the value to return 10665305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 10675305534fSZach Atkins 10685305534fSZach Atkins Level: intermediate 10695305534fSZach Atkins 10705305534fSZach Atkins Notes: 1071af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 10725305534fSZach Atkins 10735305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 10749314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`. 10755305534fSZach Atkins 10765305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 10775305534fSZach Atkins 10785305534fSZach Atkins See `PetscOptionsEList()` for when the choices are given in a string array 10795305534fSZach Atkins 10805305534fSZach Atkins To get a listing of all currently specified options, 10815305534fSZach Atkins see `PetscOptionsView()` or `PetscOptionsGetAll()` 10825305534fSZach Atkins 1083af27ebaaSBarry Smith Developer Note: 10845305534fSZach Atkins This cannot check for invalid selection because of things like `MATAIJ` that are not included in the list 10855305534fSZach Atkins 10865305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 10875305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 10885305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 10895305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 10905305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 10915305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 10925305534fSZach Atkins M*/ 10931ff8fb82SZach Atkins #define PetscOptionsFList(opt, ltext, man, list, currentvalue, value, len, set) PetscOptionsFList_Private(PetscOptionsObject, opt, ltext, man, list, currentvalue, value, len, set) 10945305534fSZach Atkins 10955305534fSZach Atkins /*MC 10965305534fSZach Atkins PetscOptionsEList - Puts a list of option values that a single one may be selected from 10975305534fSZach Atkins 10985305534fSZach Atkins Synopsis: 10995305534fSZach Atkins #include <petscoptions.h> 11005305534fSZach Atkins PetscErrorCode PetscOptionsEList(const char opt[], const char ltext[], const char man[], const char *const *list, PetscInt ntext, const char currentvalue[], PetscInt *value, PetscBool *set) 11015305534fSZach Atkins 11025305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11035305534fSZach Atkins 11045305534fSZach Atkins Input Parameters: 11055305534fSZach Atkins + opt - option name 11065305534fSZach Atkins . ltext - short string that describes the option 11075305534fSZach Atkins . man - manual page with additional information on option 11085305534fSZach Atkins . list - the possible choices (one of these must be selected, anything else is invalid) 11095305534fSZach Atkins . ntext - number of choices 11105305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 11115305534fSZach Atkins .vb 11129314d9b7SBarry Smith PetscOptionsEList(..., obj->value,&value,&set); 11139314d9b7SBarry Smith .ve if (set) { 11145305534fSZach Atkins 11155305534fSZach Atkins Output Parameters: 11165305534fSZach Atkins + value - the index of the value to return 11175305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11185305534fSZach Atkins 11195305534fSZach Atkins Level: intermediate 11205305534fSZach Atkins 11215305534fSZach Atkins Notes: 1122af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 11235305534fSZach Atkins 11245305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 11259314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`. 11265305534fSZach Atkins 11275305534fSZach Atkins See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList()` 11285305534fSZach Atkins 11295305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 11305305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 11315305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 11325305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 11335305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 11345305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEnum()` 11355305534fSZach Atkins M*/ 11361ff8fb82SZach Atkins #define PetscOptionsEList(opt, ltext, man, list, ntext, currentvalue, value, set) PetscOptionsEList_Private(PetscOptionsObject, opt, ltext, man, list, ntext, currentvalue, value, set) 11375305534fSZach Atkins 11385305534fSZach Atkins /*MC 11395305534fSZach Atkins PetscOptionsRealArray - Gets an array of double values for a particular 11405305534fSZach Atkins option in the database. The values must be separated with commas with 11415305534fSZach Atkins no intervening spaces. 11425305534fSZach Atkins 11435305534fSZach Atkins Synopsis: 11445305534fSZach Atkins #include <petscoptions.h> 11455305534fSZach Atkins PetscErrorCode PetscOptionsRealArray(const char opt[], const char text[], const char man[], PetscReal value[], PetscInt *n, PetscBool *set) 11465305534fSZach Atkins 11475305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11485305534fSZach Atkins 11495305534fSZach Atkins Input Parameters: 11505305534fSZach Atkins + opt - the option one is seeking 11515305534fSZach Atkins . text - short string describing option 11525305534fSZach Atkins . man - manual page for option 11535305534fSZach Atkins - n - maximum number of values that value has room for 11545305534fSZach Atkins 11555305534fSZach Atkins Output Parameters: 11565305534fSZach Atkins + value - location to copy values 11575305534fSZach Atkins . n - actual number of values found 11585305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11595305534fSZach Atkins 11605305534fSZach Atkins Level: beginner 11615305534fSZach Atkins 11625305534fSZach Atkins Note: 1163af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 11645305534fSZach Atkins 11655305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 11665305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 11675305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 11685305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 11695305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 11705305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 11715305534fSZach Atkins M*/ 11721ff8fb82SZach Atkins #define PetscOptionsRealArray(opt, text, man, currentvalue, value, set) PetscOptionsRealArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 11735305534fSZach Atkins 11745305534fSZach Atkins /*MC 11755305534fSZach Atkins PetscOptionsScalarArray - Gets an array of `PetscScalar` values for a particular 11765305534fSZach Atkins option in the database. The values must be separated with commas with 11775305534fSZach Atkins no intervening spaces. 11785305534fSZach Atkins 11795305534fSZach Atkins Synopsis: 11805305534fSZach Atkins #include <petscoptions.h> 11815305534fSZach Atkins PetscErrorCode PetscOptionsScalarArray(const char opt[], const char text[], const char man[], PetscScalar value[], PetscInt *n, PetscBool *set) 11825305534fSZach Atkins 11835305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11845305534fSZach Atkins 11855305534fSZach Atkins Input Parameters: 11865305534fSZach Atkins + opt - the option one is seeking 11875305534fSZach Atkins . text - short string describing option 11885305534fSZach Atkins . man - manual page for option 11895305534fSZach Atkins - n - maximum number of values allowed in the value array 11905305534fSZach Atkins 11915305534fSZach Atkins Output Parameters: 11925305534fSZach Atkins + value - location to copy values 11935305534fSZach Atkins . n - actual number of values found 11945305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11955305534fSZach Atkins 11965305534fSZach Atkins Level: beginner 11975305534fSZach Atkins 11985305534fSZach Atkins Note: 1199af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 12005305534fSZach Atkins 12015305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 12025305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 12035305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12045305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12055305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12065305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12075305534fSZach Atkins M*/ 12081ff8fb82SZach Atkins #define PetscOptionsScalarArray(opt, text, man, currentvalue, value, set) PetscOptionsScalarArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12095305534fSZach Atkins 12105305534fSZach Atkins /*MC 12115305534fSZach Atkins PetscOptionsIntArray - Gets an array of integers for a particular 12125305534fSZach Atkins option in the database. 12135305534fSZach Atkins 12145305534fSZach Atkins Synopsis: 12155305534fSZach Atkins #include <petscoptions.h> 12165305534fSZach Atkins PetscErrorCode PetscOptionsIntArray(const char opt[], const char text[], const char man[], PetscInt value[], PetscInt *n, PetscBool *set) 12175305534fSZach Atkins 12185305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 12195305534fSZach Atkins 12205305534fSZach Atkins Input Parameters: 12215305534fSZach Atkins + opt - the option one is seeking 12225305534fSZach Atkins . text - short string describing option 12235305534fSZach Atkins . man - manual page for option 12245305534fSZach Atkins - n - maximum number of values 12255305534fSZach Atkins 12265305534fSZach Atkins Output Parameters: 12275305534fSZach Atkins + value - location to copy values 12285305534fSZach Atkins . n - actual number of values found 12295305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 12305305534fSZach Atkins 12315305534fSZach Atkins Level: beginner 12325305534fSZach Atkins 12335305534fSZach Atkins Notes: 12345305534fSZach Atkins The array can be passed as 12355305534fSZach Atkins + a comma separated list - 0,1,2,3,4,5,6,7 12365305534fSZach Atkins . a range (start\-end+1) - 0-8 12375305534fSZach Atkins . a range with given increment (start\-end+1:inc) - 0-7:2 12385305534fSZach Atkins - a combination of values and ranges separated by commas - 0,1-8,8-15:2 12395305534fSZach Atkins 12405305534fSZach Atkins There must be no intervening spaces between the values. 12415305534fSZach Atkins 1242af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 12435305534fSZach Atkins 12445305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 12455305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 12465305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12475305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12485305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12495305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12505305534fSZach Atkins M*/ 12511ff8fb82SZach Atkins #define PetscOptionsIntArray(opt, text, man, currentvalue, value, set) PetscOptionsIntArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12525305534fSZach Atkins 12535305534fSZach Atkins /*MC 12545305534fSZach Atkins PetscOptionsStringArray - Gets an array of string values for a particular 12555305534fSZach Atkins option in the database. The values must be separated with commas with 12565305534fSZach Atkins no intervening spaces. 12575305534fSZach Atkins 12585305534fSZach Atkins Synopsis: 12595305534fSZach Atkins #include <petscoptions.h> 12605305534fSZach Atkins PetscErrorCode PetscOptionsStringArray(const char opt[], const char text[], const char man[], char *value[], PetscInt *nmax, PetscBool *set) 12615305534fSZach Atkins 12625305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()`; No Fortran Support 12635305534fSZach Atkins 12645305534fSZach Atkins Input Parameters: 12655305534fSZach Atkins + opt - the option one is seeking 12665305534fSZach Atkins . text - short string describing option 12675305534fSZach Atkins . man - manual page for option 12685305534fSZach Atkins - nmax - maximum number of strings 12695305534fSZach Atkins 12705305534fSZach Atkins Output Parameters: 12715305534fSZach Atkins + value - location to copy strings 12725305534fSZach Atkins . nmax - actual number of strings found 12735305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 12745305534fSZach Atkins 12755305534fSZach Atkins Level: beginner 12765305534fSZach Atkins 12775305534fSZach Atkins Notes: 12785305534fSZach Atkins The user should pass in an array of pointers to char, to hold all the 12795305534fSZach Atkins strings returned by this function. 12805305534fSZach Atkins 12815305534fSZach Atkins The user is responsible for deallocating the strings that are 12825305534fSZach Atkins returned. 12835305534fSZach Atkins 1284af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 12855305534fSZach Atkins 12865305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 12875305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 12885305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12895305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12905305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12915305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12925305534fSZach Atkins M*/ 12931ff8fb82SZach Atkins #define PetscOptionsStringArray(opt, text, man, currentvalue, value, set) PetscOptionsStringArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12945305534fSZach Atkins 12955305534fSZach Atkins /*MC 12965305534fSZach Atkins PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 12975305534fSZach Atkins option in the database. The values must be separated with commas with 12985305534fSZach Atkins no intervening spaces. 12995305534fSZach Atkins 13005305534fSZach Atkins Synopsis: 13015305534fSZach Atkins #include <petscoptions.h> 13025305534fSZach Atkins PetscErrorCode PetscOptionsBoolArray(const char opt[], const char text[], const char man[], PetscBool value[], PetscInt *n, PetscBool *set) 13035305534fSZach Atkins 13045305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 13055305534fSZach Atkins 13065305534fSZach Atkins Input Parameters: 13075305534fSZach Atkins + opt - the option one is seeking 13085305534fSZach Atkins . text - short string describing option 13095305534fSZach Atkins . man - manual page for option 13105305534fSZach Atkins - n - maximum number of values allowed in the value array 13115305534fSZach Atkins 13125305534fSZach Atkins Output Parameters: 13135305534fSZach Atkins + value - location to copy values 13145305534fSZach Atkins . n - actual number of values found 13155305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 13165305534fSZach Atkins 13175305534fSZach Atkins Level: beginner 13185305534fSZach Atkins 13195305534fSZach Atkins Notes: 13205305534fSZach Atkins The user should pass in an array of `PetscBool` 13215305534fSZach Atkins 1322af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 13235305534fSZach Atkins 13245305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 13255305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 13265305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 13275305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 13285305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 13295305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 13305305534fSZach Atkins M*/ 13311ff8fb82SZach Atkins #define PetscOptionsBoolArray(opt, text, man, currentvalue, value, set) PetscOptionsBoolArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 13325305534fSZach Atkins 13335305534fSZach Atkins /*MC 13345305534fSZach Atkins PetscOptionsEnumArray - Gets an array of enum values for a particular 13355305534fSZach Atkins option in the database. 13365305534fSZach Atkins 13375305534fSZach Atkins Synopsis: 13385305534fSZach Atkins #include <petscoptions.h> 13395305534fSZach Atkins PetscErrorCode PetscOptionsEnumArray(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum value[], PetscInt *n, PetscBool *set) 13405305534fSZach Atkins 13415305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 13425305534fSZach Atkins 13435305534fSZach Atkins Input Parameters: 13445305534fSZach Atkins + opt - the option one is seeking 13455305534fSZach Atkins . text - short string describing option 13465305534fSZach Atkins . man - manual page for option 13475305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 13485305534fSZach Atkins - n - maximum number of values allowed in the value array 13495305534fSZach Atkins 13505305534fSZach Atkins Output Parameters: 13515305534fSZach Atkins + value - location to copy values 13525305534fSZach Atkins . n - actual number of values found 13535305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 13545305534fSZach Atkins 13555305534fSZach Atkins Level: beginner 13565305534fSZach Atkins 13575305534fSZach Atkins Notes: 13585305534fSZach Atkins The array must be passed as a comma separated list. 13595305534fSZach Atkins 13605305534fSZach Atkins There must be no intervening spaces between the values. 13615305534fSZach Atkins 1362af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 13635305534fSZach Atkins 13645305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 13655305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 13665305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 13675305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 13685305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 13695305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 13705305534fSZach Atkins M*/ 13711ff8fb82SZach Atkins #define PetscOptionsEnumArray(opt, text, man, list, value, n, set) PetscOptionsEnumArray_Private(PetscOptionsObject, opt, text, man, list, value, n, set) 13725305534fSZach Atkins 13735305534fSZach Atkins /*MC 13745305534fSZach Atkins PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with `newname` 13755305534fSZach Atkins 13765305534fSZach Atkins Prints a deprecation warning, unless an option is supplied to suppress. 13775305534fSZach Atkins 13785305534fSZach Atkins Logically Collective 13795305534fSZach Atkins 13805305534fSZach Atkins Input Parameters: 13815305534fSZach Atkins + oldname - the old, deprecated option 13825305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed 13835305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9" 13845305534fSZach Atkins - info - additional information string, or `NULL`. 13855305534fSZach Atkins 13865305534fSZach Atkins Options Database Key: 13875305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings 13885305534fSZach Atkins 13895305534fSZach Atkins Level: developer 13905305534fSZach Atkins 13915305534fSZach Atkins Notes: 13925305534fSZach Atkins If `newname` is provided then the options database will automatically check the database for `oldname`. 13935305534fSZach Atkins 13945305534fSZach Atkins The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the 13955305534fSZach Atkins new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call. 13965305534fSZach Atkins See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails. 13975305534fSZach Atkins 13985305534fSZach Atkins Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`. 1399af27ebaaSBarry Smith Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information 14005305534fSZach Atkins If newname is provided, the old option is replaced. Otherwise, it remains in the options database. 14015305534fSZach Atkins If an option is not replaced, the info argument should be used to advise the user on how to proceed. 14025305534fSZach Atkins There is a limit on the length of the warning printed, so very long strings provided as info may be truncated. 14035305534fSZach Atkins 14045305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()` 14055305534fSZach Atkins M*/ 14061ff8fb82SZach Atkins #define PetscOptionsDeprecated(opt, text, man, info) PetscOptionsDeprecated_Private(PetscOptionsObject, opt, text, man, info) 14075305534fSZach Atkins 14085305534fSZach Atkins /*MC 14095305534fSZach Atkins PetscOptionsDeprecatedMoObject - mark an option as deprecated in the global PetscOptionsObject, optionally replacing it with `newname` 14105305534fSZach Atkins 14115305534fSZach Atkins Prints a deprecation warning, unless an option is supplied to suppress. 14125305534fSZach Atkins 14135305534fSZach Atkins Logically Collective 14145305534fSZach Atkins 14155305534fSZach Atkins Input Parameters: 14165305534fSZach Atkins + oldname - the old, deprecated option 14175305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed 14185305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9" 14195305534fSZach Atkins - info - additional information string, or `NULL`. 14205305534fSZach Atkins 14215305534fSZach Atkins Options Database Key: 14225305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings 14235305534fSZach Atkins 14245305534fSZach Atkins Level: developer 14255305534fSZach Atkins 14265305534fSZach Atkins Notes: 14275305534fSZach Atkins If `newname` is provided then the options database will automatically check the database for `oldname`. 14285305534fSZach Atkins 14295305534fSZach Atkins The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the 14305305534fSZach Atkins new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call. 14315305534fSZach Atkins See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails. 14325305534fSZach Atkins 1433af27ebaaSBarry Smith Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information 14345305534fSZach Atkins If newname is provided, the old option is replaced. Otherwise, it remains in the options database. 14355305534fSZach Atkins If an option is not replaced, the info argument should be used to advise the user on how to proceed. 14365305534fSZach Atkins There is a limit on the length of the warning printed, so very long strings provided as info may be truncated. 14375305534fSZach Atkins 14385305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()` 14395305534fSZach Atkins M*/ 14401ff8fb82SZach Atkins #define PetscOptionsDeprecatedNoObject(opt, text, man, info) PetscOptionsDeprecated_Private(NULL, opt, text, man, info) 14415f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1442e55864a3SBarry Smith 14434416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum, PetscEnum *, PetscBool *); 14445a856986SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt, PetscInt *, PetscBool *, PetscInt, PetscInt); 14456497c311SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsMPIInt_Private(PetscOptionItems *, const char[], const char[], const char[], PetscMPIInt, PetscMPIInt *, PetscBool *, PetscMPIInt, PetscMPIInt); 144652ce0ab5SPierre Jolivet PETSC_EXTERN PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal, PetscReal *, PetscBool *, PetscReal, PetscReal); 14474416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar, PetscScalar *, PetscBool *); 14484416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsName_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14494416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsString_Private(PetscOptionItems *, const char[], const char[], const char[], const char[], char *, size_t, PetscBool *); 14504416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool, PetscBool *, PetscBool *); 14514416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14524416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14534416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14544416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFList_Private(PetscOptionItems *, const char[], const char[], const char[], PetscFunctionList, const char[], char[], size_t, PetscBool *); 14554416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEList_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscInt, const char[], PetscInt *, PetscBool *); 14564416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal[], PetscInt *, PetscBool *); 14574416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar[], PetscInt *, PetscBool *); 14584416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt[], PetscInt *, PetscBool *); 14594416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *, const char[], const char[], const char[], char *[], PetscInt *, PetscBool *); 14604416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool[], PetscInt *, PetscBool *); 14614416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum[], PetscInt *, PetscBool *); 14629f3a6782SPatrick Sanan PETSC_EXTERN PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *, const char[], const char[], const char[], const char[]); 1463cffb1e40SBarry Smith 1464e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSAWsDestroy(void); 1465f8d0b74dSMatthew Knepley 1466dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject, PetscErrorCode (*)(PetscObject, PetscOptionItems *, void *), PetscErrorCode (*)(PetscObject, void *), void *); 1467dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject, PetscOptionItems *); 1468447ac60bSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject); 1469447ac60bSBarry Smith 1470f4bc716fSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeftError(void); 1471