xref: /petsc/include/petscoptions.h (revision 1690c2ae071c7584458d4e437df7b47bc4686b3c)
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 *);
312d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetInt(PetscOptions, const char[], const char[], PetscInt *, PetscBool *);
326497c311SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetMPIInt(PetscOptions, const char[], const char[], PetscMPIInt *, PetscBool *);
332d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnum(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscBool *);
342d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEList(PetscOptions, const char[], const char[], const char *const *, PetscInt, PetscInt *, PetscBool *);
35c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetReal(PetscOptions, const char[], const char[], PetscReal *, PetscBool *);
36c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalar(PetscOptions, const char[], const char[], PetscScalar *, PetscBool *);
372d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetString(PetscOptions, const char[], const char[], char[], size_t, PetscBool *);
382d747510SLisandro Dalcin 
392d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetBoolArray(PetscOptions, const char[], const char[], PetscBool[], PetscInt *, PetscBool *);
402d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnumArray(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscInt *, PetscBool *);
41c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetIntArray(PetscOptions, const char[], const char[], PetscInt[], PetscInt *, PetscBool *);
42c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetRealArray(PetscOptions, const char[], const char[], PetscReal[], PetscInt *, PetscBool *);
43c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalarArray(PetscOptions, const char[], const char[], PetscScalar[], PetscInt *, PetscBool *);
44c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetStringArray(PetscOptions, const char[], const char[], char *[], PetscInt *, PetscBool *);
453a3b2205SBarry Smith 
462d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsValidKey(const char[], PetscBool *);
47c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetAlias(PetscOptions, const char[], const char[]);
48c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetValue(PetscOptions, const char[], const char[]);
49c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClearValue(PetscOptions, const char[]);
502d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPair(PetscOptions, const char[], const char[], const char *[], PetscBool *);
513a3b2205SBarry Smith 
522d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetAll(PetscOptions, char *[]);
53c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsAllUsed(PetscOptions, PetscInt *);
542d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsUsed(PetscOptions, const char[], PetscBool *);
55c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeft(PetscOptions);
561ab23d95SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftGet(PetscOptions, PetscInt *, char ***, char ***);
575b191818SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftRestore(PetscOptions, PetscInt *, char ***, char ***);
58c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsView(PetscOptions, PetscViewer);
594b0e389bSBarry Smith 
602d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsReject(PetscOptions, const char[], const char[], const char[]);
61c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsert(PetscOptions, int *, char ***, const char[]);
62c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertFile(MPI_Comm, PetscOptions, const char[], PetscBool);
635c23ca1cSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertFileYAML(MPI_Comm, PetscOptions, const char[], PetscBool);
64c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertString(PetscOptions, const char[]);
65080f0011SToby Isaac PETSC_EXTERN PetscErrorCode PetscOptionsInsertStringYAML(PetscOptions, const char[]);
66d06005cbSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertArgs(PetscOptions, int, char **);
67c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClear(PetscOptions);
68c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPush(PetscOptions, const char[]);
69c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPop(PetscOptions);
705d0dffe5SBarry Smith 
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsGetenv(MPI_Comm, const char[], char[], size_t, PetscBool *);
722d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToBool(const char[], PetscBool *);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToInt(const char[], PetscInt *);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToReal(const char[], PetscReal *);
752d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToScalar(const char[], PetscScalar *);
762e8a6d31SBarry Smith 
779355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*)(const char[], const char[], PetscOptionSource, void *), void *, PetscErrorCode (*)(void **));
789355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorDefault(const char[], const char[], PetscOptionSource, void *);
79081c24baSBoyana Norris 
8084761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectSetOptions(PetscObject, PetscOptions);
8184761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectGetOptions(PetscObject, PetscOptions *);
8284761bfeSJed Brown 
83014dd563SJed Brown PETSC_EXTERN PetscBool PetscOptionsPublish;
84e55864a3SBarry Smith 
85e55864a3SBarry Smith /*
86e55864a3SBarry Smith     See manual page for PetscOptionsBegin()
874416b707SBarry Smith 
884416b707SBarry Smith     PetscOptionsItem and PetscOptionsItems are a single option (such as ksp_type) and a collection of such single
894416b707SBarry Smith   options being handled with a PetscOptionsBegin/End()
904416b707SBarry Smith 
91e55864a3SBarry Smith */
929371c9d4SSatish Balay typedef enum {
939371c9d4SSatish Balay   OPTION_INT,
949371c9d4SSatish Balay   OPTION_BOOL,
959371c9d4SSatish Balay   OPTION_REAL,
969371c9d4SSatish Balay   OPTION_FLIST,
979371c9d4SSatish Balay   OPTION_STRING,
989371c9d4SSatish Balay   OPTION_REAL_ARRAY,
999371c9d4SSatish Balay   OPTION_SCALAR_ARRAY,
1009371c9d4SSatish Balay   OPTION_HEAD,
1019371c9d4SSatish Balay   OPTION_INT_ARRAY,
1029371c9d4SSatish Balay   OPTION_ELIST,
1039371c9d4SSatish Balay   OPTION_BOOL_ARRAY,
1049371c9d4SSatish Balay   OPTION_STRING_ARRAY
1059371c9d4SSatish Balay } PetscOptionType;
1069355ec05SMatthew G. Knepley 
1074416b707SBarry Smith typedef struct _n_PetscOptionItem *PetscOptionItem;
1084416b707SBarry Smith struct _n_PetscOptionItem {
109e55864a3SBarry Smith   char              *option;
110e55864a3SBarry Smith   char              *text;
111e55864a3SBarry Smith   void              *data;  /* used to hold the default value and then any value it is changed to by GUI */
112e55864a3SBarry Smith   PetscFunctionList  flist; /* used for available values for PetscOptionsList() */
113e55864a3SBarry Smith   const char *const *list;  /* used for available values for PetscOptionsEList() */
114e55864a3SBarry Smith   char               nlist; /* number of entries in list */
115e55864a3SBarry Smith   char              *man;
1166497c311SBarry Smith   PetscInt           arraylength; /* number of entries in data in the case that it is an array (of PetscInt etc), never a giant value */
117e55864a3SBarry Smith   PetscBool          set;         /* the user has changed this value in the GUI */
118e55864a3SBarry Smith   PetscOptionType    type;
1194416b707SBarry Smith   PetscOptionItem    next;
120e55864a3SBarry Smith   char              *pman;
121e55864a3SBarry Smith   void              *edata;
122e55864a3SBarry Smith };
123e55864a3SBarry Smith 
1244416b707SBarry Smith typedef struct _p_PetscOptionItems {
125e55864a3SBarry Smith   PetscInt        count;
1264416b707SBarry Smith   PetscOptionItem next;
127e55864a3SBarry Smith   char           *prefix, *pprefix;
128e55864a3SBarry Smith   char           *title;
129e55864a3SBarry Smith   MPI_Comm        comm;
130e55864a3SBarry Smith   PetscBool       printhelp, changedmethod, alreadyprinted;
131e55864a3SBarry Smith   PetscObject     object;
132c5929fdfSBarry Smith   PetscOptions    options;
1334416b707SBarry Smith } PetscOptionItems;
13430de9b25SBarry Smith 
1355f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
13694bad497SJacob Faibussowitsch extern PetscOptionItems *PetscOptionsObject; /* declare this so that the PetscOptions stubs work */
1375f80ce2aSJacob Faibussowitsch PetscErrorCode           PetscOptionsBegin(MPI_Comm, const char *, const char *, const char *);
1385f80ce2aSJacob Faibussowitsch PetscErrorCode           PetscObjectOptionsBegin(PetscObject);
1395f80ce2aSJacob Faibussowitsch PetscErrorCode           PetscOptionsEnd(void);
1405f80ce2aSJacob Faibussowitsch #else
14130de9b25SBarry Smith   /*MC
14230de9b25SBarry Smith     PetscOptionsBegin - Begins a set of queries on the options database that are related and should be
1431957e957SBarry Smith      displayed on the same window of a GUI that allows the user to set the options interactively. Often one should
14416a05f60SBarry Smith      use `PetscObjectOptionsBegin()` rather than this call.
14530de9b25SBarry Smith 
146f2ba6396SBarry Smith     Synopsis:
147aaa7dc30SBarry Smith     #include <petscoptions.h>
148f2ba6396SBarry Smith     PetscErrorCode PetscOptionsBegin(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
14930de9b25SBarry Smith 
150fb455bf4SPatrick Sanan     Collective
15130de9b25SBarry Smith 
15230de9b25SBarry Smith     Input Parameters:
15330de9b25SBarry Smith +   comm - communicator that shares GUI
15476280437SVaclav Hapla .   prefix - options prefix for all options displayed on window (optional)
15530de9b25SBarry Smith .   title - short descriptive text, for example "Krylov Solver Options"
15687497f52SBarry Smith -   mansec - section of manual pages for options, for example `KSP` (optional)
15730de9b25SBarry Smith 
15830de9b25SBarry Smith     Level: intermediate
15930de9b25SBarry Smith 
16095452b02SPatrick Sanan     Notes:
161d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
162d0609cedSBarry Smith 
16387497f52SBarry Smith     The set of queries needs to be ended by a call to `PetscOptionsEnd()`.
164fb455bf4SPatrick Sanan 
16587497f52SBarry Smith     One can add subheadings with `PetscOptionsHeadBegin()`.
16630de9b25SBarry Smith 
16795452b02SPatrick Sanan     Developer Notes:
16816a05f60SBarry Smith     `PetscOptionsPublish` is set in `PetscOptionsCheckInitial_Private()` with `-saws_options`. When `PetscOptionsPublish` is set the
16916a05f60SBarry Smith     loop between `PetscOptionsBegin()` and `PetscOptionsEnd()` is run THREE times with `PetscOptionsPublishCount` of values -1,0,1.
17016a05f60SBarry Smith      Otherwise the loop is run ONCE with a `PetscOptionsPublishCount` of 1.
17116a05f60SBarry Smith +      \-1 - `PetscOptionsInt()` etc. just call `PetscOptionsGetInt()` etc.
17216a05f60SBarry Smith .      0   - The GUI objects are created in `PetscOptionsInt()` etc. and displayed in `PetscOptionsEnd()` and the options
17316a05f60SBarry Smith               database updated with user changes; `PetscOptionsGetInt()` etc. are also called.
17416a05f60SBarry Smith -      1   - `PetscOptionsInt()` etc. again call `PetscOptionsGetInt()` etc. (possibly getting new values), in addition the help message and
175fb455bf4SPatrick Sanan               default values are printed if -help was given.
17616a05f60SBarry Smith      When `PetscOptionsObject.changedmethod` is set this causes `PetscOptionsPublishCount` to be reset to -2 (so in the next loop iteration it is -1)
17716a05f60SBarry Smith      and the whole process is repeated. This is to handle when, for example, the `KSPType` is changed thus changing the list of
17816a05f60SBarry Smith      options available so they need to be redisplayed so the user can change the. Changing `PetscOptionsObjects.changedmethod` is never
179fb455bf4SPatrick Sanan      currently set.
180aee2cecaSBarry Smith 
1814bb2516aSBarry Smith      Fortran Note:
1824bb2516aSBarry Smith      Returns ierr error code per PETSc Fortran API
1834bb2516aSBarry Smith 
184db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
185db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
186db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
187db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
188c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
189db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
190db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()`
19130de9b25SBarry Smith M*/
1929371c9d4SSatish Balay   #define PetscOptionsBegin(comm, prefix, mess, sec) \
1939371c9d4SSatish Balay     do { \
1944416b707SBarry Smith       PetscOptionItems  PetscOptionsObjectBase; \
1954416b707SBarry Smith       PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \
196d0609cedSBarry Smith       PetscCall(PetscMemzero(PetscOptionsObject, sizeof(*PetscOptionsObject))); \
197e55864a3SBarry Smith       for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \
1989566063dSJacob Faibussowitsch         PetscCall(PetscOptionsBegin_Private(PetscOptionsObject, comm, prefix, mess, sec))
19930de9b25SBarry Smith 
2005fefd1ebSJed Brown   /*MC
2015fefd1ebSJed Brown     PetscObjectOptionsBegin - Begins a set of queries on the options database that are related and should be
2025fefd1ebSJed Brown     displayed on the same window of a GUI that allows the user to set the options interactively.
2035fefd1ebSJed Brown 
204f2ba6396SBarry Smith     Synopsis:
205aaa7dc30SBarry Smith     #include <petscoptions.h>
206f2ba6396SBarry Smith     PetscErrorCode PetscObjectOptionsBegin(PetscObject obj)
2075fefd1ebSJed Brown 
208c3339decSBarry Smith     Collective
2095fefd1ebSJed Brown 
2102fe279fdSBarry Smith     Input Parameter:
2115fefd1ebSJed Brown .   obj - object to set options for
2125fefd1ebSJed Brown 
2135fefd1ebSJed Brown     Level: intermediate
2145fefd1ebSJed Brown 
21595452b02SPatrick Sanan     Notes:
216d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
217d0609cedSBarry Smith 
21887497f52SBarry Smith     Needs to be ended by a call the `PetscOptionsEnd()`
219d0609cedSBarry Smith 
22087497f52SBarry Smith     Can add subheadings with `PetscOptionsHeadBegin()`
2215fefd1ebSJed Brown 
222db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
223db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
224db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
225db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
226c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
227db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
228db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2295fefd1ebSJed Brown M*/
2309371c9d4SSatish Balay   #define PetscObjectOptionsBegin(obj) \
2319371c9d4SSatish Balay     do { \
2324416b707SBarry Smith       PetscOptionItems  PetscOptionsObjectBase; \
2334416b707SBarry Smith       PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \
234c5929fdfSBarry Smith       PetscOptionsObject->options          = ((PetscObject)obj)->options; \
235e55864a3SBarry Smith       for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \
236dbbe0bcdSBarry Smith         PetscCall(PetscObjectOptionsBegin_Private(obj, PetscOptionsObject))
2373194b578SJed Brown 
23830de9b25SBarry Smith   /*MC
23930de9b25SBarry Smith     PetscOptionsEnd - Ends a set of queries on the options database that are related and should be
24030de9b25SBarry Smith      displayed on the same window of a GUI that allows the user to set the options interactively.
24130de9b25SBarry Smith 
242f2ba6396SBarry Smith     Synopsis:
243aaa7dc30SBarry Smith      #include <petscoptions.h>
244f2ba6396SBarry Smith      PetscErrorCode PetscOptionsEnd(void)
24530de9b25SBarry Smith 
2467cdbe19fSJose E. Roman     Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()`
2477cdbe19fSJose E. Roman 
24830de9b25SBarry Smith     Level: intermediate
24930de9b25SBarry Smith 
25095452b02SPatrick Sanan     Notes:
25187497f52SBarry Smith     Needs to be preceded by a call to `PetscOptionsBegin()` or `PetscObjectOptionsBegin()`
25230de9b25SBarry Smith 
253d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
254d0609cedSBarry Smith 
2554bb2516aSBarry Smith     Fortran Note:
2564bb2516aSBarry Smith     Returns ierr error code per PETSc Fortran API
2574bb2516aSBarry Smith 
258db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
259db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
260db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
261db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsHeadBegin()`,
262c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
263db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
264db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()`
26530de9b25SBarry Smith M*/
2669371c9d4SSatish Balay   #define PetscOptionsEnd() \
2679371c9d4SSatish Balay     PetscCall(PetscOptionsEnd_Private(PetscOptionsObject)); \
2689371c9d4SSatish Balay     } \
2699371c9d4SSatish Balay     } \
2709371c9d4SSatish Balay     while (0)
2715f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
27230de9b25SBarry Smith 
2734416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *, MPI_Comm, const char[], const char[], const char[]);
274dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject, PetscOptionItems *);
2754416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *);
276d0609cedSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHeadBegin(PetscOptionItems *, const char[]);
27730de9b25SBarry Smith 
2785f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
2799371c9d4SSatish Balay template <typename... T>
2809371c9d4SSatish Balay void PetscOptionsHeadBegin(T...);
281d0609cedSBarry Smith void PetscOptionsHeadEnd(void);
2829371c9d4SSatish Balay template <typename... T>
2839371c9d4SSatish Balay PetscErrorCode PetscOptionsEnum(T...);
2849371c9d4SSatish Balay template <typename... T>
2859371c9d4SSatish Balay PetscErrorCode PetscOptionsInt(T...);
2869371c9d4SSatish Balay template <typename... T>
2879371c9d4SSatish Balay PetscErrorCode PetscOptionsBoundedInt(T...);
2889371c9d4SSatish Balay template <typename... T>
2899371c9d4SSatish Balay PetscErrorCode PetscOptionsRangeInt(T...);
2909371c9d4SSatish Balay template <typename... T>
2919371c9d4SSatish Balay PetscErrorCode PetscOptionsReal(T...);
2929371c9d4SSatish Balay template <typename... T>
2939371c9d4SSatish Balay PetscErrorCode PetscOptionsScalar(T...);
2949371c9d4SSatish Balay template <typename... T>
2959371c9d4SSatish Balay PetscErrorCode PetscOptionsName(T...);
2969371c9d4SSatish Balay template <typename... T>
2979371c9d4SSatish Balay PetscErrorCode PetscOptionsString(T...);
2989371c9d4SSatish Balay template <typename... T>
2999371c9d4SSatish Balay PetscErrorCode PetscOptionsBool(T...);
3009371c9d4SSatish Balay template <typename... T>
3019371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(T...);
3029371c9d4SSatish Balay template <typename... T>
3039371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroup(T...);
3049371c9d4SSatish Balay template <typename... T>
3059371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(T...);
3069371c9d4SSatish Balay template <typename... T>
3079371c9d4SSatish Balay PetscErrorCode PetscOptionsFList(T...);
3089371c9d4SSatish Balay template <typename... T>
3099371c9d4SSatish Balay PetscErrorCode PetscOptionsEList(T...);
3109371c9d4SSatish Balay template <typename... T>
3119371c9d4SSatish Balay PetscErrorCode PetscOptionsRealArray(T...);
3129371c9d4SSatish Balay template <typename... T>
3139371c9d4SSatish Balay PetscErrorCode PetscOptionsScalarArray(T...);
3149371c9d4SSatish Balay template <typename... T>
3159371c9d4SSatish Balay PetscErrorCode PetscOptionsIntArray(T...);
3169371c9d4SSatish Balay template <typename... T>
3179371c9d4SSatish Balay PetscErrorCode PetscOptionsStringArray(T...);
3189371c9d4SSatish Balay template <typename... T>
3199371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolArray(T...);
3209371c9d4SSatish Balay template <typename... T>
3219371c9d4SSatish Balay PetscErrorCode PetscOptionsEnumArray(T...);
3229371c9d4SSatish Balay template <typename... T>
3239371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecated(T...);
3249371c9d4SSatish Balay template <typename... T>
3259371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecatedNoObject(T...);
3265f80ce2aSJacob Faibussowitsch #else
32730de9b25SBarry Smith   /*MC
328d0609cedSBarry Smith    PetscOptionsHeadBegin - Puts a heading before listing any more published options. Used, for example,
32916a05f60SBarry Smith    in `KSPSetFromOptions_GMRES()`.
330d0609cedSBarry Smith 
33116a05f60SBarry Smith    Logically Collective on the communicator passed in `PetscOptionsBegin()`
332d0609cedSBarry Smith 
333d0609cedSBarry Smith    Input Parameter:
334d0609cedSBarry Smith .  head - the heading text
335d0609cedSBarry Smith 
33687497f52SBarry Smith    Level: developer
337d0609cedSBarry Smith 
338d0609cedSBarry Smith    Notes:
339d0609cedSBarry Smith    Handles errors directly, hence does not return an error code
340d0609cedSBarry Smith 
34116a05f60SBarry Smith    Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`, and `PetscOptionsObject` created in `PetscOptionsBegin()` should be the first argument
342d0609cedSBarry Smith 
34387497f52SBarry Smith    Must be followed by a call to `PetscOptionsHeadEnd()` in the same function.
344d0609cedSBarry Smith 
345db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
346db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
347db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
348c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
349db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
350db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
351a54bb2a9SBarry Smith M*/
3529371c9d4SSatish Balay   #define PetscOptionsHeadBegin(PetscOptionsObject, head) \
3539371c9d4SSatish Balay     do { \
35448a46eb9SPierre Jolivet       if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, "  %s\n", head)); \
355d0609cedSBarry Smith     } while (0)
356d0609cedSBarry Smith 
357edd03b47SJacob Faibussowitsch   #define PetscOptionsHead(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadBegin()", ) PetscOptionsHeadBegin(__VA_ARGS__)
358d0609cedSBarry Smith 
359d0609cedSBarry Smith   /*MC
36087497f52SBarry Smith      PetscOptionsHeadEnd - Ends a section of options begun with `PetscOptionsHeadBegin()`
36116a05f60SBarry Smith      See, for example, `KSPSetFromOptions_GMRES()`.
36230de9b25SBarry Smith 
363f2ba6396SBarry Smith      Synopsis:
364aaa7dc30SBarry Smith      #include <petscoptions.h>
365d0609cedSBarry Smith      PetscErrorCode PetscOptionsHeadEnd(void)
36630de9b25SBarry Smith 
3677cdbe19fSJose E. Roman      Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()`
3687cdbe19fSJose E. Roman 
36930de9b25SBarry Smith      Level: intermediate
37030de9b25SBarry Smith 
37195452b02SPatrick Sanan      Notes:
37287497f52SBarry Smith      Must be between a `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` and a `PetscOptionsEnd()`
37330de9b25SBarry Smith 
37487497f52SBarry Smith      Must be preceded by a call to `PetscOptionsHeadBegin()` in the same function.
37530de9b25SBarry Smith 
37687497f52SBarry Smith      This needs to be used only if the code below `PetscOptionsHeadEnd()` can be run ONLY once.
37716a05f60SBarry Smith      See, for example, `PCSetFromOptions_Composite()`. This is a `return(0)` in it for early exit
378b52f573bSBarry Smith      from the function.
379b52f573bSBarry Smith 
38056752e42SBarry Smith      This is only for use with the PETSc options GUI
381b52f573bSBarry Smith 
382db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
383db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
384db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
385c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
386db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
387db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()`
38830de9b25SBarry Smith M*/
3899371c9d4SSatish Balay   #define PetscOptionsHeadEnd() \
3909371c9d4SSatish Balay     do { \
3913ba16761SJacob Faibussowitsch       if (PetscOptionsObject->count != 1) PetscFunctionReturn(PETSC_SUCCESS); \
3929371c9d4SSatish Balay     } while (0)
393d0609cedSBarry Smith 
394edd03b47SJacob Faibussowitsch   #define PetscOptionsTail(...)                                                     PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadEnd()", ) PetscOptionsHeadEnd(__VA_ARGS__)
395186905e3SBarry Smith 
3965305534fSZach Atkins /*MC
3975305534fSZach Atkins   PetscOptionsEnum - Gets the enum value for a particular option in the database.
3985305534fSZach Atkins 
3995305534fSZach Atkins   Synopsis:
4005305534fSZach Atkins   #include <petscoptions.h>
4015305534fSZach Atkins   PetscErrorCode PetscOptionsEnum(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum currentvalue, PetscEnum *value, PetscBool *set)
4025305534fSZach Atkins 
4035305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4045305534fSZach Atkins 
4055305534fSZach Atkins   Input Parameters:
4065305534fSZach Atkins + opt          - option name
4075305534fSZach Atkins . text         - short string that describes the option
4085305534fSZach Atkins . man          - manual page with additional information on option
4095305534fSZach Atkins . list         - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
4105305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
4115305534fSZach Atkins .vb
4125305534fSZach Atkins                  PetscOptionsEnum(..., obj->value,&object->value,...) or
4135305534fSZach Atkins                  value = defaultvalue
4149314d9b7SBarry Smith                  PetscOptionsEnum(..., value,&value,&set);
4159314d9b7SBarry Smith                  if (set) {
4165305534fSZach Atkins .ve
4175305534fSZach Atkins 
4185305534fSZach Atkins   Output Parameters:
4195305534fSZach Atkins + value - the  value to return
4205305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
4215305534fSZach Atkins 
4225305534fSZach Atkins   Level: beginner
4235305534fSZach Atkins 
4245305534fSZach Atkins   Notes:
4255305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
4265305534fSZach Atkins 
4279314d9b7SBarry Smith   `list` is usually something like `PCASMTypes` or some other predefined list of enum names
4285305534fSZach Atkins 
4295305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
4309314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
4315305534fSZach Atkins 
4325305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
4335305534fSZach Atkins 
4345305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
4355305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
4365305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
4375305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
4385305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
4395305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
4405305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
4415305534fSZach Atkins M*/
4421ff8fb82SZach Atkins   #define PetscOptionsEnum(opt, text, man, list, currentvalue, value, set)          PetscOptionsEnum_Private(PetscOptionsObject, opt, text, man, list, currentvalue, value, set)
4435305534fSZach Atkins 
4445305534fSZach Atkins /*MC
4455305534fSZach Atkins   PetscOptionsInt - Gets the integer value for a particular option in the database.
4465305534fSZach Atkins 
4475305534fSZach Atkins   Synopsis:
4485305534fSZach Atkins   #include <petscoptions.h>
4495305534fSZach Atkins   PetscErrorCode PetscOptionsInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set)
4505305534fSZach Atkins 
4515305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4525305534fSZach Atkins 
4535305534fSZach Atkins   Input Parameters:
4545305534fSZach Atkins + opt          - option name
4555305534fSZach Atkins . text         - short string that describes the option
4565305534fSZach Atkins . man          - manual page with additional information on option
4575305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
4585305534fSZach Atkins .vb
4595305534fSZach Atkins                  PetscOptionsInt(..., obj->value, &obj->value, ...) or
4605305534fSZach Atkins                  value = defaultvalue
4619314d9b7SBarry Smith                  PetscOptionsInt(..., value, &value, &set);
4629314d9b7SBarry Smith                  if (set) {
4635305534fSZach Atkins .ve
4645305534fSZach Atkins 
4655305534fSZach Atkins   Output Parameters:
4665305534fSZach Atkins + value - the integer value to return
4675305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
4685305534fSZach Atkins 
4695305534fSZach Atkins   Level: beginner
4705305534fSZach Atkins 
4715305534fSZach Atkins   Notes:
4725305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
4739314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
4745305534fSZach Atkins 
4755305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
4765305534fSZach Atkins 
4775305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
4785305534fSZach Atkins 
4795305534fSZach Atkins .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
4805305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
4815305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
4825305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
4835305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
4845305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
48552ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
4865305534fSZach Atkins M*/
487*1690c2aeSBarry Smith   #define PetscOptionsInt(opt, text, man, currentvalue, value, set)                 PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_INT_MIN, PETSC_INT_MAX)
4885305534fSZach Atkins 
4895305534fSZach Atkins /*MC
4906497c311SBarry Smith   PetscOptionsMPIInt - Gets the MPI integer value for a particular option in the database.
4916497c311SBarry Smith 
4926497c311SBarry Smith   Synopsis:
4936497c311SBarry Smith   #include <petscoptions.h>
4946497c311SBarry Smith   PetscErrorCode PetscOptionsMPIInt(const char opt[], const char text[], const char man[], PetscMPIInt currentvalue, PetscMPIInt *value, PetscBool *set)
4956497c311SBarry Smith 
4966497c311SBarry Smith   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4976497c311SBarry Smith 
4986497c311SBarry Smith   Input Parameters:
4996497c311SBarry Smith + opt          - option name
5006497c311SBarry Smith . text         - short string that describes the option
5016497c311SBarry Smith . man          - manual page with additional information on option
5026497c311SBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5036497c311SBarry Smith .vb
5046497c311SBarry Smith                  PetscOptionsInt(..., obj->value, &obj->value, ...) or
5056497c311SBarry Smith                  value = defaultvalue
5066497c311SBarry Smith                  PetscOptionsInt(..., value, &value, &set);
5076497c311SBarry Smith                  if (set) {
5086497c311SBarry Smith .ve
5096497c311SBarry Smith 
5106497c311SBarry Smith   Output Parameters:
5116497c311SBarry Smith + value - the MPI integer value to return
5126497c311SBarry Smith - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
5136497c311SBarry Smith 
5146497c311SBarry Smith   Level: beginner
5156497c311SBarry Smith 
5166497c311SBarry Smith   Notes:
5176497c311SBarry Smith   If the user does not supply the option at all `value` is NOT changed. Thus
5186497c311SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
5196497c311SBarry Smith 
5206497c311SBarry Smith   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
5216497c311SBarry Smith 
5226497c311SBarry Smith   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
5236497c311SBarry Smith 
5246497c311SBarry Smith .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
5256497c311SBarry Smith           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
5266497c311SBarry Smith           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
5276497c311SBarry Smith           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
5286497c311SBarry Smith           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
5296497c311SBarry Smith           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
5306497c311SBarry Smith           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
5316497c311SBarry Smith M*/
5326497c311SBarry 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)
5336497c311SBarry Smith 
5346497c311SBarry Smith /*MC
53552ce0ab5SPierre Jolivet    PetscOptionsBoundedInt - Gets an integer value greater than or equal to a given bound for a particular option in the database.
5365305534fSZach Atkins 
5375305534fSZach Atkins    Synopsis:
5385305534fSZach Atkins    #include <petscoptions.h>
5399314d9b7SBarry Smith    PetscErrorCode  PetscOptionsBoundedInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt bound)
5405305534fSZach Atkins 
5415305534fSZach Atkins    Logically Collective on the communicator passed in `PetscOptionsBegin()`
5425305534fSZach Atkins 
5435305534fSZach Atkins    Input Parameters:
5445305534fSZach Atkins +  opt          - option name
5455305534fSZach Atkins .  text         - short string that describes the option
5465305534fSZach Atkins .  man          - manual page with additional information on option
5475305534fSZach Atkins .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5485305534fSZach Atkins .vb
54952ce0ab5SPierre Jolivet   PetscOptionsBoundedInt(..., obj->value, &obj->value, ...)
5505305534fSZach Atkins .ve
5515305534fSZach Atkins or
5525305534fSZach Atkins .vb
5535305534fSZach Atkins   value = defaultvalue
55452ce0ab5SPierre Jolivet   PetscOptionsBoundedInt(..., value, &value, &set, ...);
5559314d9b7SBarry Smith   if (set) {
5565305534fSZach Atkins .ve
55752ce0ab5SPierre Jolivet -  bound - the requested value should be greater than or equal to this bound or an error is generated
5585305534fSZach Atkins 
5595305534fSZach Atkins    Output Parameters:
5605305534fSZach Atkins +  value - the integer value to return
5619314d9b7SBarry Smith -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
5625305534fSZach Atkins 
5635305534fSZach Atkins    Level: beginner
5645305534fSZach Atkins 
5655305534fSZach Atkins    Notes:
5665305534fSZach Atkins    If the user does not supply the option at all `value` is NOT changed. Thus
5679314d9b7SBarry Smith    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
5685305534fSZach Atkins 
5695305534fSZach Atkins    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
5705305534fSZach Atkins 
571af27ebaaSBarry Smith    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
5725305534fSZach Atkins 
5735305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
5745305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
5755305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
5765305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
5775305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
5785305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
57952ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
5805305534fSZach Atkins M*/
581*1690c2aeSBarry Smith   #define PetscOptionsBoundedInt(opt, text, man, currentvalue, value, set, lb)      PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_INT_MAX)
5825305534fSZach Atkins 
5835305534fSZach Atkins /*MC
5845305534fSZach Atkins    PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database.
5855305534fSZach Atkins 
5865305534fSZach Atkins    Synopsis:
5875305534fSZach Atkins    #include <petscoptions.h>
5889314d9b7SBarry Smith    PetscErrorCode PetscOptionsRangeInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt lb, PetscInt ub)
5895305534fSZach Atkins 
5905305534fSZach Atkins    Logically Collective on the communicator passed in `PetscOptionsBegin()`
5915305534fSZach Atkins 
5925305534fSZach Atkins    Input Parameters:
5935305534fSZach Atkins +  opt          - option name
5945305534fSZach Atkins .  text         - short string that describes the option
5955305534fSZach Atkins .  man          - manual page with additional information on option
5965305534fSZach Atkins .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5975305534fSZach Atkins .vb
59852ce0ab5SPierre Jolivet   PetscOptionsRangeInt(..., obj->value, &obj->value, ...)
59952ce0ab5SPierre Jolivet .ve
60052ce0ab5SPierre Jolivet or
60152ce0ab5SPierre Jolivet .vb
6025305534fSZach Atkins   value = defaultvalue
60352ce0ab5SPierre Jolivet   PetscOptionsRangeInt(..., value, &value, &set, ...);
6049314d9b7SBarry Smith   if (set) {
6055305534fSZach Atkins .ve
6065305534fSZach Atkins .  lb - the lower bound, provided value must be greater than or equal to this value or an error is generated
6075305534fSZach Atkins -  ub - the upper bound, provided value must be less than or equal to this value or an error is generated
6085305534fSZach Atkins 
6095305534fSZach Atkins    Output Parameters:
6105305534fSZach Atkins +  value - the integer value to return
6119314d9b7SBarry Smith -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
6125305534fSZach Atkins 
6135305534fSZach Atkins    Level: beginner
6145305534fSZach Atkins 
6155305534fSZach Atkins    Notes:
6165305534fSZach Atkins    If the user does not supply the option at all `value` is NOT changed. Thus
6179314d9b7SBarry Smith    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
6185305534fSZach Atkins 
6195305534fSZach Atkins    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
6205305534fSZach Atkins 
621af27ebaaSBarry Smith    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
6225305534fSZach Atkins 
6235305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
6245305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()`
6255305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
6265305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
6275305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
6285305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
62952ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
6305305534fSZach Atkins M*/
6311ff8fb82SZach Atkins   #define PetscOptionsRangeInt(opt, text, man, currentvalue, value, set, lb, ub)    PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub)
6325305534fSZach Atkins 
6335305534fSZach Atkins /*MC
63452ce0ab5SPierre Jolivet   PetscOptionsReal - Gets a `PetscReal` value for a particular option in the database.
6355305534fSZach Atkins 
6365305534fSZach Atkins   Synopsis:
6375305534fSZach Atkins   #include <petscoptions.h>
6385305534fSZach Atkins   PetscErrorCode PetscOptionsReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set)
6395305534fSZach Atkins 
6405305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
6415305534fSZach Atkins 
6425305534fSZach Atkins   Input Parameters:
6435305534fSZach Atkins + opt          - option name
6445305534fSZach Atkins . text         - short string that describes the option
6455305534fSZach Atkins . man          - manual page with additional information on option
6465305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6475305534fSZach Atkins .vb
6485305534fSZach Atkins                  PetscOptionsReal(..., obj->value,&obj->value,...) or
6495305534fSZach Atkins                  value = defaultvalue
6509314d9b7SBarry Smith                  PetscOptionsReal(..., value,&value,&set);
6519314d9b7SBarry Smith                  if (set) {
6525305534fSZach Atkins .ve
6535305534fSZach Atkins 
6545305534fSZach Atkins   Output Parameters:
6555305534fSZach Atkins + value - the value to return
6565305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
6575305534fSZach Atkins 
6585305534fSZach Atkins   Level: beginner
6595305534fSZach Atkins 
6605305534fSZach Atkins   Notes:
6615305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
6629314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
6635305534fSZach Atkins 
6645305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
6655305534fSZach Atkins 
666af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
6675305534fSZach Atkins 
6685305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
6695305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
6705305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
6715305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
6725305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
6735305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67452ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
6755305534fSZach Atkins M*/
67652ce0ab5SPierre Jolivet   #define PetscOptionsReal(opt, text, man, currentvalue, value, set)                PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MIN_REAL, PETSC_MAX_REAL)
67752ce0ab5SPierre Jolivet 
67852ce0ab5SPierre Jolivet /*MC
67952ce0ab5SPierre Jolivet    PetscOptionsBoundedReal - Gets a `PetscReal` value greater than or equal to a given bound for a particular option in the database.
68052ce0ab5SPierre Jolivet 
68152ce0ab5SPierre Jolivet    Synopsis:
68252ce0ab5SPierre Jolivet    #include <petscoptions.h>
68352ce0ab5SPierre Jolivet    PetscErrorCode  PetscOptionsBoundedReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal bound)
68452ce0ab5SPierre Jolivet 
68552ce0ab5SPierre Jolivet    Logically Collective on the communicator passed in `PetscOptionsBegin()`
68652ce0ab5SPierre Jolivet 
68752ce0ab5SPierre Jolivet    Input Parameters:
68852ce0ab5SPierre Jolivet +  opt          - option name
68952ce0ab5SPierre Jolivet .  text         - short string that describes the option
69052ce0ab5SPierre Jolivet .  man          - manual page with additional information on option
69152ce0ab5SPierre Jolivet .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
69252ce0ab5SPierre Jolivet .vb
69352ce0ab5SPierre Jolivet   PetscOptionsBoundedReal(..., obj->value, &obj->value, ...)
69452ce0ab5SPierre Jolivet .ve
69552ce0ab5SPierre Jolivet or
69652ce0ab5SPierre Jolivet .vb
69752ce0ab5SPierre Jolivet   value = defaultvalue
69852ce0ab5SPierre Jolivet   PetscOptionsBoundedReal(..., value, &value, &set, ...);
69952ce0ab5SPierre Jolivet   if (set) {
70052ce0ab5SPierre Jolivet .ve
70152ce0ab5SPierre Jolivet -  bound - the requested value should be greater than or equal to this bound or an error is generated
70252ce0ab5SPierre Jolivet 
70352ce0ab5SPierre Jolivet    Output Parameters:
70452ce0ab5SPierre Jolivet +  value - the real value to return
70552ce0ab5SPierre Jolivet -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
70652ce0ab5SPierre Jolivet 
70752ce0ab5SPierre Jolivet    Level: beginner
70852ce0ab5SPierre Jolivet 
70952ce0ab5SPierre Jolivet    Notes:
71052ce0ab5SPierre Jolivet    If the user does not supply the option at all `value` is NOT changed. Thus
71152ce0ab5SPierre Jolivet    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
71252ce0ab5SPierre Jolivet 
71352ce0ab5SPierre Jolivet    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
71452ce0ab5SPierre Jolivet 
71552ce0ab5SPierre Jolivet    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
71652ce0ab5SPierre Jolivet 
71752ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
71852ce0ab5SPierre Jolivet           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
71952ce0ab5SPierre Jolivet           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
72052ce0ab5SPierre Jolivet           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
72152ce0ab5SPierre Jolivet           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
72252ce0ab5SPierre Jolivet           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
72352ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedInt()`, `PetscOptionsRangeReal()`
72452ce0ab5SPierre Jolivet M*/
72552ce0ab5SPierre Jolivet   #define PetscOptionsBoundedReal(opt, text, man, currentvalue, value, set, lb)     PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_MAX_REAL)
72652ce0ab5SPierre Jolivet 
72752ce0ab5SPierre Jolivet /*MC
72852ce0ab5SPierre Jolivet    PetscOptionsRangeReal - Gets a `PetscReal` value within a range of values for a particular option in the database.
72952ce0ab5SPierre Jolivet 
73052ce0ab5SPierre Jolivet    Synopsis:
73152ce0ab5SPierre Jolivet    #include <petscoptions.h>
73252ce0ab5SPierre Jolivet    PetscErrorCode PetscOptionsRangeReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal lb, PetscReal ub)
73352ce0ab5SPierre Jolivet 
73452ce0ab5SPierre Jolivet    Logically Collective on the communicator passed in `PetscOptionsBegin()`
73552ce0ab5SPierre Jolivet 
73652ce0ab5SPierre Jolivet    Input Parameters:
73752ce0ab5SPierre Jolivet +  opt          - option name
73852ce0ab5SPierre Jolivet .  text         - short string that describes the option
73952ce0ab5SPierre Jolivet .  man          - manual page with additional information on option
74052ce0ab5SPierre Jolivet .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
74152ce0ab5SPierre Jolivet .vb
74252ce0ab5SPierre Jolivet   PetscOptionsRangeReal(..., obj->value, &obj->value, ...)
74352ce0ab5SPierre Jolivet .ve
74452ce0ab5SPierre Jolivet or
74552ce0ab5SPierre Jolivet .vb
74652ce0ab5SPierre Jolivet   value = defaultvalue
74752ce0ab5SPierre Jolivet   PetscOptionsRangeReal(..., value, &value, &set, ...);
74852ce0ab5SPierre Jolivet   if (set) {
74952ce0ab5SPierre Jolivet .ve
75052ce0ab5SPierre Jolivet .  lb - the lower bound, provided value must be greater than or equal to this value or an error is generated
75152ce0ab5SPierre Jolivet -  ub - the upper bound, provided value must be less than or equal to this value or an error is generated
75252ce0ab5SPierre Jolivet 
75352ce0ab5SPierre Jolivet    Output Parameters:
75452ce0ab5SPierre Jolivet +  value - the value to return
75552ce0ab5SPierre Jolivet -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
75652ce0ab5SPierre Jolivet 
75752ce0ab5SPierre Jolivet    Level: beginner
75852ce0ab5SPierre Jolivet 
75952ce0ab5SPierre Jolivet    Notes:
76052ce0ab5SPierre Jolivet    If the user does not supply the option at all `value` is NOT changed. Thus
76152ce0ab5SPierre Jolivet    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
76252ce0ab5SPierre Jolivet 
76352ce0ab5SPierre Jolivet    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
76452ce0ab5SPierre Jolivet 
76552ce0ab5SPierre Jolivet    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
76652ce0ab5SPierre Jolivet 
76752ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
76852ce0ab5SPierre Jolivet           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()`
76952ce0ab5SPierre Jolivet           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
77052ce0ab5SPierre Jolivet           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
77152ce0ab5SPierre Jolivet           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
77252ce0ab5SPierre Jolivet           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
77352ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRangeInt()`, `PetscOptionsBoundedReal()`
77452ce0ab5SPierre Jolivet M*/
77552ce0ab5SPierre Jolivet   #define PetscOptionsRangeReal(opt, text, man, currentvalue, value, set, lb, ub)   PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub)
7765305534fSZach Atkins 
7775305534fSZach Atkins /*MC
7785305534fSZach Atkins   PetscOptionsScalar - Gets the `PetscScalar` value for a particular option in the database.
7795305534fSZach Atkins 
7805305534fSZach Atkins   Synopsis:
7815305534fSZach Atkins   #include <petscoptions.h>
7825305534fSZach Atkins   PetscErrorCode PetscOptionsScalar(const char opt[], const char text[], const char man[], PetscScalar currentvalue, PetscScalar *value, PetscBool *set)
7835305534fSZach Atkins 
7845305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
7855305534fSZach Atkins 
7865305534fSZach Atkins   Input Parameters:
7875305534fSZach Atkins + opt          - option name
7885305534fSZach Atkins . text         - short string that describes the option
7895305534fSZach Atkins . man          - manual page with additional information on option
7905305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7915305534fSZach Atkins .vb
7925305534fSZach Atkins                  PetscOptionsScalar(..., obj->value,&obj->value,...) or
7935305534fSZach Atkins                  value = defaultvalue
7949314d9b7SBarry Smith                  PetscOptionsScalar(..., value,&value,&set);
7959314d9b7SBarry Smith                  if (set) {
7965305534fSZach Atkins .ve
7975305534fSZach Atkins 
7985305534fSZach Atkins   Output Parameters:
7995305534fSZach Atkins + value - the value to return
8005305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
8015305534fSZach Atkins 
8025305534fSZach Atkins   Level: beginner
8035305534fSZach Atkins 
8045305534fSZach Atkins   Notes:
8055305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
8069314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
8075305534fSZach Atkins 
8085305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
8095305534fSZach Atkins 
810af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8115305534fSZach Atkins 
8125305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8135305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8145305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8155305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8165305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8175305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8185305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8195305534fSZach Atkins M*/
8201ff8fb82SZach Atkins   #define PetscOptionsScalar(opt, text, man, currentvalue, value, set)              PetscOptionsScalar_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
8215305534fSZach Atkins 
8225305534fSZach Atkins /*MC
8235305534fSZach 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
8245305534fSZach Atkins   its value is set to false.
8255305534fSZach Atkins 
8265305534fSZach Atkins   Synopsis:
8275305534fSZach Atkins   #include <petscoptions.h>
8289314d9b7SBarry Smith   PetscErrorCode PetscOptionsName(const char opt[], const char text[], const char man[], PetscBool *set)
8295305534fSZach Atkins 
8305305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
8315305534fSZach Atkins 
8325305534fSZach Atkins   Input Parameters:
8335305534fSZach Atkins + opt  - option name
8345305534fSZach Atkins . text - short string that describes the option
8355305534fSZach Atkins - man  - manual page with additional information on option
8365305534fSZach Atkins 
8375305534fSZach Atkins   Output Parameter:
8389314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
8395305534fSZach Atkins 
8405305534fSZach Atkins   Level: beginner
8415305534fSZach Atkins 
8425305534fSZach Atkins   Note:
843af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8445305534fSZach Atkins 
8455305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8465305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8475305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8485305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8495305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8505305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8515305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8525305534fSZach Atkins M*/
8539314d9b7SBarry Smith   #define PetscOptionsName(opt, text, man, set)                                     PetscOptionsName_Private(PetscOptionsObject, opt, text, man, set)
8545305534fSZach Atkins 
8555305534fSZach Atkins /*MC
8565305534fSZach Atkins   PetscOptionsString - Gets the string value for a particular option in the database.
8575305534fSZach Atkins 
8585305534fSZach Atkins   Synopsis:
8595305534fSZach Atkins   #include <petscoptions.h>
8605305534fSZach Atkins   PetscErrorCode PetscOptionsString(const char opt[], const char text[], const char man[], const char currentvalue[], char value[], size_t len, PetscBool *set)
8615305534fSZach Atkins 
8625305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
8635305534fSZach Atkins 
8645305534fSZach Atkins   Input Parameters:
8655305534fSZach Atkins + opt          - option name
8665305534fSZach Atkins . text         - short string that describes the option
8675305534fSZach Atkins . man          - manual page with additional information on option
8685305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
8695305534fSZach Atkins - len          - length of the result string including null terminator
8705305534fSZach Atkins 
8715305534fSZach Atkins   Output Parameters:
8725305534fSZach Atkins + value - the value to return
8735305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
8745305534fSZach Atkins 
8755305534fSZach Atkins   Level: beginner
8765305534fSZach Atkins 
8775305534fSZach Atkins   Notes:
878af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8795305534fSZach Atkins 
8809314d9b7SBarry Smith   If the user provided no string (for example `-optionname` `-someotheroption`) `set` is set to `PETSC_TRUE` (and the string is filled with nulls).
8815305534fSZach Atkins 
8825305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
8839314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
8845305534fSZach Atkins 
8855305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
8865305534fSZach Atkins 
8875305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8885305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8895305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8905305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8915305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8925305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8935305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8945305534fSZach Atkins M*/
8951ff8fb82SZach Atkins   #define PetscOptionsString(opt, text, man, currentvalue, value, len, set)         PetscOptionsString_Private(PetscOptionsObject, opt, text, man, currentvalue, value, len, set)
8965305534fSZach Atkins 
8975305534fSZach Atkins /*MC
8985305534fSZach Atkins   PetscOptionsBool - Determines if a particular option is in the database with a true or false
8995305534fSZach Atkins 
9005305534fSZach Atkins   Synopsis:
9015305534fSZach Atkins   #include <petscoptions.h>
9025305534fSZach Atkins   PetscErrorCode PetscOptionsBool(const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool *flg, PetscBool *set)
9035305534fSZach Atkins 
9045305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
9055305534fSZach Atkins 
9065305534fSZach Atkins   Input Parameters:
9075305534fSZach Atkins + opt          - option name
9085305534fSZach Atkins . text         - short string that describes the option
9095305534fSZach Atkins . man          - manual page with additional information on option
9105305534fSZach Atkins - currentvalue - the current value
9115305534fSZach Atkins 
9125305534fSZach Atkins   Output Parameters:
9135305534fSZach Atkins + flg - `PETSC_TRUE` or `PETSC_FALSE`
9145305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
9155305534fSZach Atkins 
9165305534fSZach Atkins   Level: beginner
9175305534fSZach Atkins 
9185305534fSZach Atkins   Notes:
9195305534fSZach Atkins   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
9205305534fSZach Atkins   FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
9215305534fSZach Atkins 
9229314d9b7SBarry 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`
9235305534fSZach Atkins   is equivalent to `-requested_bool true`
9245305534fSZach Atkins 
9255305534fSZach Atkins   If the user does not supply the option at all `flg` is NOT changed. Thus
9269314d9b7SBarry Smith   you should ALWAYS initialize the `flg` variable if you access it without first checking that the `set` flag is `PETSC_TRUE`.
9275305534fSZach Atkins 
9285305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
9295305534fSZach Atkins 
9305305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
9315305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
9325305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
9335305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
9345305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
9355305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
9365305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
9375305534fSZach Atkins M*/
9381ff8fb82SZach Atkins   #define PetscOptionsBool(opt, text, man, currentvalue, value, set)                PetscOptionsBool_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
9395305534fSZach Atkins 
9405305534fSZach Atkins /*MC
9415305534fSZach Atkins   PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
9425305534fSZach Atkins   which at most a single value can be true.
9435305534fSZach Atkins 
9445305534fSZach Atkins   Synopsis:
9455305534fSZach Atkins   #include <petscoptions.h>
9469314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[], const char text[], const char man[], PetscBool *set)
9475305534fSZach Atkins 
9485305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
9495305534fSZach Atkins 
9505305534fSZach Atkins   Input Parameters:
9515305534fSZach Atkins + opt  - option name
9525305534fSZach Atkins . text - short string that describes the option
9535305534fSZach Atkins - man  - manual page with additional information on option
9545305534fSZach Atkins 
9555305534fSZach Atkins   Output Parameter:
9569314d9b7SBarry Smith . set - whether that option was set or not
9575305534fSZach Atkins 
9585305534fSZach Atkins   Level: intermediate
9595305534fSZach Atkins 
9605305534fSZach Atkins   Notes:
961af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
9625305534fSZach Atkins 
9635305534fSZach Atkins   Must be followed by 0 or more `PetscOptionsBoolGroup()`s and `PetscOptionsBoolGroupEnd()`
9645305534fSZach Atkins 
9655305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
9665305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
9675305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
9685305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
9695305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
9705305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
9715305534fSZach Atkins M*/
9729314d9b7SBarry Smith   #define PetscOptionsBoolGroupBegin(opt, text, man, set)                           PetscOptionsBoolGroupBegin_Private(PetscOptionsObject, opt, text, man, set)
9735305534fSZach Atkins 
9745305534fSZach Atkins /*MC
9755305534fSZach Atkins   PetscOptionsBoolGroup - One in a series of logical queries on the options database for
9765305534fSZach Atkins   which at most a single value can be true.
9775305534fSZach Atkins 
9785305534fSZach Atkins   Synopsis:
9795305534fSZach Atkins   #include <petscoptions.h>
9809314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroup(const char opt[], const char text[], const char man[], PetscBool *set)
9815305534fSZach Atkins 
9825305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
9835305534fSZach Atkins 
9845305534fSZach Atkins   Input Parameters:
9855305534fSZach Atkins + opt  - option name
9865305534fSZach Atkins . text - short string that describes the option
9875305534fSZach Atkins - man  - manual page with additional information on option
9885305534fSZach Atkins 
9895305534fSZach Atkins   Output Parameter:
9909314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
9915305534fSZach Atkins 
9925305534fSZach Atkins   Level: intermediate
9935305534fSZach Atkins 
9945305534fSZach Atkins   Notes:
995af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
9965305534fSZach Atkins 
9975305534fSZach Atkins   Must follow a `PetscOptionsBoolGroupBegin()` and preceded a `PetscOptionsBoolGroupEnd()`
9985305534fSZach Atkins 
9995305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10005305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10015305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10025305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10035305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10045305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
10055305534fSZach Atkins M*/
10069314d9b7SBarry Smith   #define PetscOptionsBoolGroup(opt, text, man, set)                                PetscOptionsBoolGroup_Private(PetscOptionsObject, opt, text, man, set)
10075305534fSZach Atkins 
10085305534fSZach Atkins /*MC
10095305534fSZach Atkins   PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
10105305534fSZach Atkins   which at most a single value can be true.
10115305534fSZach Atkins 
10125305534fSZach Atkins   Synopsis:
10135305534fSZach Atkins   #include <petscoptions.h>
10149314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[], const char text[], const char man[], PetscBool  *set)
10155305534fSZach Atkins 
10165305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
10175305534fSZach Atkins 
10185305534fSZach Atkins   Input Parameters:
10195305534fSZach Atkins + opt  - option name
10205305534fSZach Atkins . text - short string that describes the option
10215305534fSZach Atkins - man  - manual page with additional information on option
10225305534fSZach Atkins 
10235305534fSZach Atkins   Output Parameter:
10249314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
10255305534fSZach Atkins 
10265305534fSZach Atkins   Level: intermediate
10275305534fSZach Atkins 
10285305534fSZach Atkins   Notes:
1029af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
10305305534fSZach Atkins 
10315305534fSZach Atkins   Must follow a `PetscOptionsBoolGroupBegin()`
10325305534fSZach Atkins 
10335305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10345305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10355305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10365305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10375305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10385305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
10395305534fSZach Atkins M*/
10409314d9b7SBarry Smith   #define PetscOptionsBoolGroupEnd(opt, text, man, set)                             PetscOptionsBoolGroupEnd_Private(PetscOptionsObject, opt, text, man, set)
10415305534fSZach Atkins 
10425305534fSZach Atkins /*MC
10435305534fSZach Atkins   PetscOptionsFList - Puts a list of option values that a single one may be selected from
10445305534fSZach Atkins 
10455305534fSZach Atkins   Synopsis:
10465305534fSZach Atkins   #include <petscoptions.h>
10475305534fSZach Atkins   PetscErrorCode PetscOptionsFList(const char opt[], const char ltext[], const char man[], PetscFunctionList list, const char currentvalue[], char value[], size_t len, PetscBool *set)
10485305534fSZach Atkins 
10495305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
10505305534fSZach Atkins 
10515305534fSZach Atkins   Input Parameters:
10525305534fSZach Atkins + opt          - option name
10535305534fSZach Atkins . ltext        - short string that describes the option
10545305534fSZach Atkins . man          - manual page with additional information on option
10555305534fSZach Atkins . list         - the possible choices
10565305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10575305534fSZach Atkins .vb
10589314d9b7SBarry Smith                  PetscOptionsFlist(..., obj->value,value,len,&set);
10599314d9b7SBarry Smith                  if (set) {
10605305534fSZach Atkins .ve
10615305534fSZach Atkins - len          - the length of the character array value
10625305534fSZach Atkins 
10635305534fSZach Atkins   Output Parameters:
10645305534fSZach Atkins + value - the value to return
10655305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
10665305534fSZach Atkins 
10675305534fSZach Atkins   Level: intermediate
10685305534fSZach Atkins 
10695305534fSZach Atkins   Notes:
1070af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
10715305534fSZach Atkins 
10725305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
10739314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`.
10745305534fSZach Atkins 
10755305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
10765305534fSZach Atkins 
10775305534fSZach Atkins   See `PetscOptionsEList()` for when the choices are given in a string array
10785305534fSZach Atkins 
10795305534fSZach Atkins   To get a listing of all currently specified options,
10805305534fSZach Atkins   see `PetscOptionsView()` or `PetscOptionsGetAll()`
10815305534fSZach Atkins 
1082af27ebaaSBarry Smith   Developer Note:
10835305534fSZach Atkins   This cannot check for invalid selection because of things like `MATAIJ` that are not included in the list
10845305534fSZach Atkins 
10855305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10865305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10875305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10885305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10895305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10905305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()`
10915305534fSZach Atkins M*/
10921ff8fb82SZach Atkins   #define PetscOptionsFList(opt, ltext, man, list, currentvalue, value, len, set)   PetscOptionsFList_Private(PetscOptionsObject, opt, ltext, man, list, currentvalue, value, len, set)
10935305534fSZach Atkins 
10945305534fSZach Atkins /*MC
10955305534fSZach Atkins   PetscOptionsEList - Puts a list of option values that a single one may be selected from
10965305534fSZach Atkins 
10975305534fSZach Atkins   Synopsis:
10985305534fSZach Atkins   #include <petscoptions.h>
10995305534fSZach 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)
11005305534fSZach Atkins 
11015305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
11025305534fSZach Atkins 
11035305534fSZach Atkins   Input Parameters:
11045305534fSZach Atkins + opt          - option name
11055305534fSZach Atkins . ltext        - short string that describes the option
11065305534fSZach Atkins . man          - manual page with additional information on option
11075305534fSZach Atkins . list         - the possible choices (one of these must be selected, anything else is invalid)
11085305534fSZach Atkins . ntext        - number of choices
11095305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11105305534fSZach Atkins .vb
11119314d9b7SBarry Smith                  PetscOptionsEList(..., obj->value,&value,&set);
11129314d9b7SBarry Smith .ve                 if (set) {
11135305534fSZach Atkins 
11145305534fSZach Atkins   Output Parameters:
11155305534fSZach Atkins + value - the index of the value to return
11165305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
11175305534fSZach Atkins 
11185305534fSZach Atkins   Level: intermediate
11195305534fSZach Atkins 
11205305534fSZach Atkins   Notes:
1121af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
11225305534fSZach Atkins 
11235305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
11249314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`.
11255305534fSZach Atkins 
11265305534fSZach Atkins   See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList()`
11275305534fSZach Atkins 
11285305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
11295305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
11305305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
11315305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
11325305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
11335305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEnum()`
11345305534fSZach Atkins M*/
11351ff8fb82SZach Atkins   #define PetscOptionsEList(opt, ltext, man, list, ntext, currentvalue, value, set) PetscOptionsEList_Private(PetscOptionsObject, opt, ltext, man, list, ntext, currentvalue, value, set)
11365305534fSZach Atkins 
11375305534fSZach Atkins /*MC
11385305534fSZach Atkins   PetscOptionsRealArray - Gets an array of double values for a particular
11395305534fSZach Atkins   option in the database. The values must be separated with commas with
11405305534fSZach Atkins   no intervening spaces.
11415305534fSZach Atkins 
11425305534fSZach Atkins   Synopsis:
11435305534fSZach Atkins   #include <petscoptions.h>
11445305534fSZach Atkins   PetscErrorCode PetscOptionsRealArray(const char opt[], const char text[], const char man[], PetscReal value[], PetscInt *n, PetscBool *set)
11455305534fSZach Atkins 
11465305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
11475305534fSZach Atkins 
11485305534fSZach Atkins   Input Parameters:
11495305534fSZach Atkins + opt  - the option one is seeking
11505305534fSZach Atkins . text - short string describing option
11515305534fSZach Atkins . man  - manual page for option
11525305534fSZach Atkins - n    - maximum number of values that value has room for
11535305534fSZach Atkins 
11545305534fSZach Atkins   Output Parameters:
11555305534fSZach Atkins + value - location to copy values
11565305534fSZach Atkins . n     - actual number of values found
11575305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
11585305534fSZach Atkins 
11595305534fSZach Atkins   Level: beginner
11605305534fSZach Atkins 
11615305534fSZach Atkins   Note:
1162af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
11635305534fSZach Atkins 
11645305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
11655305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
11665305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
11675305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
11685305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
11695305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
11705305534fSZach Atkins M*/
11711ff8fb82SZach Atkins   #define PetscOptionsRealArray(opt, text, man, currentvalue, value, set)           PetscOptionsRealArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
11725305534fSZach Atkins 
11735305534fSZach Atkins /*MC
11745305534fSZach Atkins   PetscOptionsScalarArray - Gets an array of `PetscScalar` values for a particular
11755305534fSZach Atkins   option in the database. The values must be separated with commas with
11765305534fSZach Atkins   no intervening spaces.
11775305534fSZach Atkins 
11785305534fSZach Atkins   Synopsis:
11795305534fSZach Atkins   #include <petscoptions.h>
11805305534fSZach Atkins   PetscErrorCode PetscOptionsScalarArray(const char opt[], const char text[], const char man[], PetscScalar value[], PetscInt *n, PetscBool *set)
11815305534fSZach Atkins 
11825305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
11835305534fSZach Atkins 
11845305534fSZach Atkins   Input Parameters:
11855305534fSZach Atkins + opt  - the option one is seeking
11865305534fSZach Atkins . text - short string describing option
11875305534fSZach Atkins . man  - manual page for option
11885305534fSZach Atkins - n    - maximum number of values allowed in the value array
11895305534fSZach Atkins 
11905305534fSZach Atkins   Output Parameters:
11915305534fSZach Atkins + value - location to copy values
11925305534fSZach Atkins . n     - actual number of values found
11935305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
11945305534fSZach Atkins 
11955305534fSZach Atkins   Level: beginner
11965305534fSZach Atkins 
11975305534fSZach Atkins   Note:
1198af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
11995305534fSZach Atkins 
12005305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12015305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12025305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12035305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12045305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12055305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12065305534fSZach Atkins M*/
12071ff8fb82SZach Atkins   #define PetscOptionsScalarArray(opt, text, man, currentvalue, value, set)         PetscOptionsScalarArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
12085305534fSZach Atkins 
12095305534fSZach Atkins /*MC
12105305534fSZach Atkins   PetscOptionsIntArray - Gets an array of integers for a particular
12115305534fSZach Atkins   option in the database.
12125305534fSZach Atkins 
12135305534fSZach Atkins   Synopsis:
12145305534fSZach Atkins   #include <petscoptions.h>
12155305534fSZach Atkins   PetscErrorCode PetscOptionsIntArray(const char opt[], const char text[], const char man[], PetscInt value[], PetscInt *n, PetscBool *set)
12165305534fSZach Atkins 
12175305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
12185305534fSZach Atkins 
12195305534fSZach Atkins   Input Parameters:
12205305534fSZach Atkins + opt  - the option one is seeking
12215305534fSZach Atkins . text - short string describing option
12225305534fSZach Atkins . man  - manual page for option
12235305534fSZach Atkins - n    - maximum number of values
12245305534fSZach Atkins 
12255305534fSZach Atkins   Output Parameters:
12265305534fSZach Atkins + value - location to copy values
12275305534fSZach Atkins . n     - actual number of values found
12285305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
12295305534fSZach Atkins 
12305305534fSZach Atkins   Level: beginner
12315305534fSZach Atkins 
12325305534fSZach Atkins   Notes:
12335305534fSZach Atkins   The array can be passed as
12345305534fSZach Atkins +   a comma separated list -                                  0,1,2,3,4,5,6,7
12355305534fSZach Atkins .   a range (start\-end+1) -                                  0-8
12365305534fSZach Atkins .   a range with given increment (start\-end+1:inc) -         0-7:2
12375305534fSZach Atkins -   a combination of values and ranges separated by commas -  0,1-8,8-15:2
12385305534fSZach Atkins 
12395305534fSZach Atkins   There must be no intervening spaces between the values.
12405305534fSZach Atkins 
1241af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
12425305534fSZach Atkins 
12435305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12445305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12455305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12465305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12475305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12485305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12495305534fSZach Atkins M*/
12501ff8fb82SZach Atkins   #define PetscOptionsIntArray(opt, text, man, currentvalue, value, set)            PetscOptionsIntArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
12515305534fSZach Atkins 
12525305534fSZach Atkins /*MC
12535305534fSZach Atkins   PetscOptionsStringArray - Gets an array of string values for a particular
12545305534fSZach Atkins   option in the database. The values must be separated with commas with
12555305534fSZach Atkins   no intervening spaces.
12565305534fSZach Atkins 
12575305534fSZach Atkins   Synopsis:
12585305534fSZach Atkins   #include <petscoptions.h>
12595305534fSZach Atkins   PetscErrorCode PetscOptionsStringArray(const char opt[], const char text[], const char man[], char *value[], PetscInt *nmax, PetscBool  *set)
12605305534fSZach Atkins 
12615305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`; No Fortran Support
12625305534fSZach Atkins 
12635305534fSZach Atkins   Input Parameters:
12645305534fSZach Atkins + opt  - the option one is seeking
12655305534fSZach Atkins . text - short string describing option
12665305534fSZach Atkins . man  - manual page for option
12675305534fSZach Atkins - nmax - maximum number of strings
12685305534fSZach Atkins 
12695305534fSZach Atkins   Output Parameters:
12705305534fSZach Atkins + value - location to copy strings
12715305534fSZach Atkins . nmax  - actual number of strings found
12725305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
12735305534fSZach Atkins 
12745305534fSZach Atkins   Level: beginner
12755305534fSZach Atkins 
12765305534fSZach Atkins   Notes:
12775305534fSZach Atkins   The user should pass in an array of pointers to char, to hold all the
12785305534fSZach Atkins   strings returned by this function.
12795305534fSZach Atkins 
12805305534fSZach Atkins   The user is responsible for deallocating the strings that are
12815305534fSZach Atkins   returned.
12825305534fSZach Atkins 
1283af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
12845305534fSZach Atkins 
12855305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12865305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12875305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12885305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12895305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12905305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12915305534fSZach Atkins M*/
12921ff8fb82SZach Atkins   #define PetscOptionsStringArray(opt, text, man, currentvalue, value, set)         PetscOptionsStringArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
12935305534fSZach Atkins 
12945305534fSZach Atkins /*MC
12955305534fSZach Atkins   PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
12965305534fSZach Atkins   option in the database. The values must be separated with commas with
12975305534fSZach Atkins   no intervening spaces.
12985305534fSZach Atkins 
12995305534fSZach Atkins   Synopsis:
13005305534fSZach Atkins   #include <petscoptions.h>
13015305534fSZach Atkins   PetscErrorCode PetscOptionsBoolArray(const char opt[], const char text[], const char man[], PetscBool value[], PetscInt *n, PetscBool *set)
13025305534fSZach Atkins 
13035305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
13045305534fSZach Atkins 
13055305534fSZach Atkins   Input Parameters:
13065305534fSZach Atkins + opt  - the option one is seeking
13075305534fSZach Atkins . text - short string describing option
13085305534fSZach Atkins . man  - manual page for option
13095305534fSZach Atkins - n    - maximum number of values allowed in the value array
13105305534fSZach Atkins 
13115305534fSZach Atkins   Output Parameters:
13125305534fSZach Atkins + value - location to copy values
13135305534fSZach Atkins . n     - actual number of values found
13145305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
13155305534fSZach Atkins 
13165305534fSZach Atkins   Level: beginner
13175305534fSZach Atkins 
13185305534fSZach Atkins   Notes:
13195305534fSZach Atkins   The user should pass in an array of `PetscBool`
13205305534fSZach Atkins 
1321af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
13225305534fSZach Atkins 
13235305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
13245305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
13255305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
13265305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
13275305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
13285305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
13295305534fSZach Atkins M*/
13301ff8fb82SZach Atkins   #define PetscOptionsBoolArray(opt, text, man, currentvalue, value, set)           PetscOptionsBoolArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
13315305534fSZach Atkins 
13325305534fSZach Atkins /*MC
13335305534fSZach Atkins   PetscOptionsEnumArray - Gets an array of enum values for a particular
13345305534fSZach Atkins   option in the database.
13355305534fSZach Atkins 
13365305534fSZach Atkins   Synopsis:
13375305534fSZach Atkins   #include <petscoptions.h>
13385305534fSZach Atkins   PetscErrorCode PetscOptionsEnumArray(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum value[], PetscInt *n, PetscBool *set)
13395305534fSZach Atkins 
13405305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
13415305534fSZach Atkins 
13425305534fSZach Atkins   Input Parameters:
13435305534fSZach Atkins + opt  - the option one is seeking
13445305534fSZach Atkins . text - short string describing option
13455305534fSZach Atkins . man  - manual page for option
13465305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
13475305534fSZach Atkins - n    - maximum number of values allowed in the value array
13485305534fSZach Atkins 
13495305534fSZach Atkins   Output Parameters:
13505305534fSZach Atkins + value - location to copy values
13515305534fSZach Atkins . n     - actual number of values found
13525305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
13535305534fSZach Atkins 
13545305534fSZach Atkins   Level: beginner
13555305534fSZach Atkins 
13565305534fSZach Atkins   Notes:
13575305534fSZach Atkins   The array must be passed as a comma separated list.
13585305534fSZach Atkins 
13595305534fSZach Atkins   There must be no intervening spaces between the values.
13605305534fSZach Atkins 
1361af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
13625305534fSZach Atkins 
13635305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
13645305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
13655305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
13665305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
13675305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
13685305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
13695305534fSZach Atkins M*/
13701ff8fb82SZach Atkins   #define PetscOptionsEnumArray(opt, text, man, list, value, n, set)                PetscOptionsEnumArray_Private(PetscOptionsObject, opt, text, man, list, value, n, set)
13715305534fSZach Atkins 
13725305534fSZach Atkins /*MC
13735305534fSZach Atkins   PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with `newname`
13745305534fSZach Atkins 
13755305534fSZach Atkins   Prints a deprecation warning, unless an option is supplied to suppress.
13765305534fSZach Atkins 
13775305534fSZach Atkins   Logically Collective
13785305534fSZach Atkins 
13795305534fSZach Atkins   Input Parameters:
13805305534fSZach Atkins + oldname - the old, deprecated option
13815305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed
13825305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9"
13835305534fSZach Atkins - info    - additional information string, or `NULL`.
13845305534fSZach Atkins 
13855305534fSZach Atkins   Options Database Key:
13865305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings
13875305534fSZach Atkins 
13885305534fSZach Atkins   Level: developer
13895305534fSZach Atkins 
13905305534fSZach Atkins   Notes:
13915305534fSZach Atkins   If `newname` is provided then the options database will automatically check the database for `oldname`.
13925305534fSZach Atkins 
13935305534fSZach Atkins   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
13945305534fSZach Atkins   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
13955305534fSZach Atkins   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
13965305534fSZach Atkins 
13975305534fSZach Atkins   Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
1398af27ebaaSBarry Smith   Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information
13995305534fSZach Atkins   If newname is provided, the old option is replaced. Otherwise, it remains in the options database.
14005305534fSZach Atkins   If an option is not replaced, the info argument should be used to advise the user on how to proceed.
14015305534fSZach Atkins   There is a limit on the length of the warning printed, so very long strings provided as info may be truncated.
14025305534fSZach Atkins 
14035305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
14045305534fSZach Atkins M*/
14051ff8fb82SZach Atkins   #define PetscOptionsDeprecated(opt, text, man, info)                              PetscOptionsDeprecated_Private(PetscOptionsObject, opt, text, man, info)
14065305534fSZach Atkins 
14075305534fSZach Atkins /*MC
14085305534fSZach Atkins   PetscOptionsDeprecatedMoObject - mark an option as deprecated in the global PetscOptionsObject, optionally replacing it with `newname`
14095305534fSZach Atkins 
14105305534fSZach Atkins   Prints a deprecation warning, unless an option is supplied to suppress.
14115305534fSZach Atkins 
14125305534fSZach Atkins   Logically Collective
14135305534fSZach Atkins 
14145305534fSZach Atkins   Input Parameters:
14155305534fSZach Atkins + oldname - the old, deprecated option
14165305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed
14175305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9"
14185305534fSZach Atkins - info    - additional information string, or `NULL`.
14195305534fSZach Atkins 
14205305534fSZach Atkins   Options Database Key:
14215305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings
14225305534fSZach Atkins 
14235305534fSZach Atkins   Level: developer
14245305534fSZach Atkins 
14255305534fSZach Atkins   Notes:
14265305534fSZach Atkins   If `newname` is provided then the options database will automatically check the database for `oldname`.
14275305534fSZach Atkins 
14285305534fSZach Atkins   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
14295305534fSZach Atkins   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
14305305534fSZach Atkins   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
14315305534fSZach Atkins 
1432af27ebaaSBarry Smith   Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information
14335305534fSZach Atkins   If newname is provided, the old option is replaced. Otherwise, it remains in the options database.
14345305534fSZach Atkins   If an option is not replaced, the info argument should be used to advise the user on how to proceed.
14355305534fSZach Atkins   There is a limit on the length of the warning printed, so very long strings provided as info may be truncated.
14365305534fSZach Atkins 
14375305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
14385305534fSZach Atkins M*/
14391ff8fb82SZach Atkins   #define PetscOptionsDeprecatedNoObject(opt, text, man, info)                      PetscOptionsDeprecated_Private(NULL, opt, text, man, info)
14405f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1441e55864a3SBarry Smith 
14424416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum, PetscEnum *, PetscBool *);
14435a856986SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt, PetscInt *, PetscBool *, PetscInt, PetscInt);
14446497c311SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsMPIInt_Private(PetscOptionItems *, const char[], const char[], const char[], PetscMPIInt, PetscMPIInt *, PetscBool *, PetscMPIInt, PetscMPIInt);
144552ce0ab5SPierre Jolivet PETSC_EXTERN PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal, PetscReal *, PetscBool *, PetscReal, PetscReal);
14464416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar, PetscScalar *, PetscBool *);
14474416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsName_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *);
14484416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsString_Private(PetscOptionItems *, const char[], const char[], const char[], const char[], char *, size_t, PetscBool *);
14494416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool, PetscBool *, PetscBool *);
14504416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *);
14514416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *);
14524416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *);
14534416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFList_Private(PetscOptionItems *, const char[], const char[], const char[], PetscFunctionList, const char[], char[], size_t, PetscBool *);
14544416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEList_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscInt, const char[], PetscInt *, PetscBool *);
14554416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal[], PetscInt *, PetscBool *);
14564416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar[], PetscInt *, PetscBool *);
14574416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt[], PetscInt *, PetscBool *);
14584416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *, const char[], const char[], const char[], char *[], PetscInt *, PetscBool *);
14594416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool[], PetscInt *, PetscBool *);
14604416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum[], PetscInt *, PetscBool *);
14619f3a6782SPatrick Sanan PETSC_EXTERN PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *, const char[], const char[], const char[], const char[]);
1462cffb1e40SBarry Smith 
1463e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSAWsDestroy(void);
1464f8d0b74dSMatthew Knepley 
1465dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject, PetscErrorCode (*)(PetscObject, PetscOptionItems *, void *), PetscErrorCode (*)(PetscObject, void *), void *);
1466dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject, PetscOptionItems *);
1467447ac60bSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1468447ac60bSBarry Smith 
1469f4bc716fSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeftError(void);
1470