xref: /petsc/include/petscoptions.h (revision fe7aa59fd2aa140db60e6d3bcaddfca3cbec1354)
1f26ada1bSBarry Smith /*
237f753daSBarry Smith    Routines to determine options set in the options database.
3f26ada1bSBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
5ac09b921SBarry Smith 
62c8e378dSBarry Smith #include <petscsys.h>
7c619b03eSJed Brown #include <petscviewertypes.h>
83a3b2205SBarry Smith 
9ac09b921SBarry Smith /* SUBMANSEC = Sys */
10ac09b921SBarry Smith 
119355ec05SMatthew G. Knepley typedef enum {
129355ec05SMatthew G. Knepley   PETSC_OPT_CODE,
139355ec05SMatthew G. Knepley   PETSC_OPT_COMMAND_LINE,
149355ec05SMatthew G. Knepley   PETSC_OPT_FILE,
159355ec05SMatthew G. Knepley   PETSC_OPT_ENVIRONMENT,
169355ec05SMatthew G. Knepley   NUM_PETSC_OPT_SOURCE
179355ec05SMatthew G. Knepley } PetscOptionSource;
189355ec05SMatthew G. Knepley 
19c5c1f447SLisandro Dalcin #define PETSC_MAX_OPTION_NAME 512
20c5929fdfSBarry Smith typedef struct _n_PetscOptions *PetscOptions;
21c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode     PetscOptionsCreate(PetscOptions *);
22b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode     PetscOptionsPush(PetscOptions);
23b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode     PetscOptionsPop(void);
24c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode     PetscOptionsDestroy(PetscOptions *);
252d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode     PetscOptionsCreateDefault(void);
262d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode     PetscOptionsDestroyDefault(void);
27c5929fdfSBarry Smith 
282d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsHasHelp(PetscOptions, PetscBool *);
29c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHasName(PetscOptions, const char[], const char[], PetscBool *);
30c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetBool(PetscOptions, const char[], const char[], PetscBool *, PetscBool *);
319e296098SJunchao Zhang PETSC_EXTERN PetscErrorCode PetscOptionsGetBool3(PetscOptions, const char[], const char[], PetscBool3 *, PetscBool *);
322d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetInt(PetscOptions, const char[], const char[], PetscInt *, PetscBool *);
336497c311SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetMPIInt(PetscOptions, const char[], const char[], PetscMPIInt *, PetscBool *);
342d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnum(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscBool *);
352d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEList(PetscOptions, const char[], const char[], const char *const *, PetscInt, PetscInt *, PetscBool *);
36c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetReal(PetscOptions, const char[], const char[], PetscReal *, PetscBool *);
37c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalar(PetscOptions, const char[], const char[], PetscScalar *, PetscBool *);
382d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetString(PetscOptions, const char[], const char[], char[], size_t, PetscBool *);
392d747510SLisandro Dalcin 
402d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetBoolArray(PetscOptions, const char[], const char[], PetscBool[], PetscInt *, PetscBool *);
412d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnumArray(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscInt *, PetscBool *);
42c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetIntArray(PetscOptions, const char[], const char[], PetscInt[], PetscInt *, PetscBool *);
43c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetRealArray(PetscOptions, const char[], const char[], PetscReal[], PetscInt *, PetscBool *);
44c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalarArray(PetscOptions, const char[], const char[], PetscScalar[], PetscInt *, PetscBool *);
45c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetStringArray(PetscOptions, const char[], const char[], char *[], PetscInt *, PetscBool *);
463a3b2205SBarry Smith 
472d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsValidKey(const char[], PetscBool *);
48c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetAlias(PetscOptions, const char[], const char[]);
49c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetValue(PetscOptions, const char[], const char[]);
50c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClearValue(PetscOptions, const char[]);
512d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPair(PetscOptions, const char[], const char[], const char *[], PetscBool *);
523a3b2205SBarry Smith 
532d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetAll(PetscOptions, char *[]);
54c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsAllUsed(PetscOptions, PetscInt *);
552d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsUsed(PetscOptions, const char[], PetscBool *);
56c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeft(PetscOptions);
571ab23d95SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftGet(PetscOptions, PetscInt *, char ***, char ***);
585b191818SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftRestore(PetscOptions, PetscInt *, char ***, char ***);
59c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsView(PetscOptions, PetscViewer);
604b0e389bSBarry Smith 
612d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsReject(PetscOptions, const char[], const char[], const char[]);
62c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsert(PetscOptions, int *, char ***, const char[]);
63c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertFile(MPI_Comm, PetscOptions, const char[], PetscBool);
645c23ca1cSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertFileYAML(MPI_Comm, PetscOptions, const char[], PetscBool);
65c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertString(PetscOptions, const char[]);
66080f0011SToby Isaac PETSC_EXTERN PetscErrorCode PetscOptionsInsertStringYAML(PetscOptions, const char[]);
67ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertArgs(PetscOptions, int, const char *const *);
68c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClear(PetscOptions);
69c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPush(PetscOptions, const char[]);
70c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPop(PetscOptions);
715d0dffe5SBarry Smith 
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsGetenv(MPI_Comm, const char[], char[], size_t, PetscBool *);
732d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToBool(const char[], PetscBool *);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToInt(const char[], PetscInt *);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToReal(const char[], PetscReal *);
762d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToScalar(const char[], PetscScalar *);
772e8a6d31SBarry Smith 
7849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*)(const char[], const char[], PetscOptionSource, void *), void *, PetscCtxDestroyFn *);
799355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorDefault(const char[], const char[], PetscOptionSource, void *);
80081c24baSBoyana Norris 
8184761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectSetOptions(PetscObject, PetscOptions);
8284761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectGetOptions(PetscObject, PetscOptions *);
8384761bfeSJed Brown 
84014dd563SJed Brown PETSC_EXTERN PetscBool PetscOptionsPublish;
85e55864a3SBarry Smith 
86e55864a3SBarry Smith /*
87e55864a3SBarry Smith     See manual page for PetscOptionsBegin()
884416b707SBarry Smith 
894416b707SBarry Smith     PetscOptionsItem and PetscOptionsItems are a single option (such as ksp_type) and a collection of such single
904416b707SBarry Smith   options being handled with a PetscOptionsBegin/End()
914416b707SBarry Smith 
92e55864a3SBarry Smith */
939371c9d4SSatish Balay typedef enum {
949371c9d4SSatish Balay   OPTION_INT,
959371c9d4SSatish Balay   OPTION_BOOL,
969371c9d4SSatish Balay   OPTION_REAL,
979371c9d4SSatish Balay   OPTION_FLIST,
989371c9d4SSatish Balay   OPTION_STRING,
999371c9d4SSatish Balay   OPTION_REAL_ARRAY,
1009371c9d4SSatish Balay   OPTION_SCALAR_ARRAY,
1019371c9d4SSatish Balay   OPTION_HEAD,
1029371c9d4SSatish Balay   OPTION_INT_ARRAY,
1039371c9d4SSatish Balay   OPTION_ELIST,
1049371c9d4SSatish Balay   OPTION_BOOL_ARRAY,
1059371c9d4SSatish Balay   OPTION_STRING_ARRAY
1069371c9d4SSatish Balay } PetscOptionType;
1079355ec05SMatthew G. Knepley 
1084416b707SBarry Smith typedef struct _n_PetscOptionItem *PetscOptionItem;
1094416b707SBarry Smith struct _n_PetscOptionItem {
110e55864a3SBarry Smith   char              *option;
111e55864a3SBarry Smith   char              *text;
112e55864a3SBarry Smith   void              *data;  /* used to hold the default value and then any value it is changed to by GUI */
113e55864a3SBarry Smith   PetscFunctionList  flist; /* used for available values for PetscOptionsList() */
114e55864a3SBarry Smith   const char *const *list;  /* used for available values for PetscOptionsEList() */
115e55864a3SBarry Smith   char               nlist; /* number of entries in list */
116e55864a3SBarry Smith   char              *man;
1176497c311SBarry Smith   PetscInt           arraylength; /* number of entries in data in the case that it is an array (of PetscInt etc), never a giant value */
118e55864a3SBarry Smith   PetscBool          set;         /* the user has changed this value in the GUI */
119e55864a3SBarry Smith   PetscOptionType    type;
1204416b707SBarry Smith   PetscOptionItem    next;
121e55864a3SBarry Smith   char              *pman;
122e55864a3SBarry Smith   void              *edata;
123e55864a3SBarry Smith };
124e55864a3SBarry Smith 
125ce78bad3SBarry Smith typedef struct _n_PetscOptionItems *PetscOptionItems;
126ce78bad3SBarry Smith struct _n_PetscOptionItems {
127e55864a3SBarry Smith   PetscInt        count;
1284416b707SBarry Smith   PetscOptionItem next;
129e55864a3SBarry Smith   char           *prefix, *pprefix;
130e55864a3SBarry Smith   char           *title;
131e55864a3SBarry Smith   MPI_Comm        comm;
132e55864a3SBarry Smith   PetscBool       printhelp, changedmethod, alreadyprinted;
133e55864a3SBarry Smith   PetscObject     object;
134c5929fdfSBarry Smith   PetscOptions    options;
135ce78bad3SBarry Smith };
13630de9b25SBarry Smith 
1375f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
138ce78bad3SBarry Smith extern PetscOptionItems PetscOptionsObject; /* declare this so that the PetscOptions stubs work */
1395f80ce2aSJacob Faibussowitsch PetscErrorCode          PetscOptionsBegin(MPI_Comm, const char *, const char *, const char *);
1405f80ce2aSJacob Faibussowitsch PetscErrorCode          PetscObjectOptionsBegin(PetscObject);
1415f80ce2aSJacob Faibussowitsch PetscErrorCode          PetscOptionsEnd(void);
1425f80ce2aSJacob Faibussowitsch #else
14330de9b25SBarry Smith   /*MC
14430de9b25SBarry Smith     PetscOptionsBegin - Begins a set of queries on the options database that are related and should be
1451957e957SBarry Smith      displayed on the same window of a GUI that allows the user to set the options interactively. Often one should
14616a05f60SBarry Smith      use `PetscObjectOptionsBegin()` rather than this call.
14730de9b25SBarry Smith 
148f2ba6396SBarry Smith     Synopsis:
149aaa7dc30SBarry Smith     #include <petscoptions.h>
150f2ba6396SBarry Smith     PetscErrorCode PetscOptionsBegin(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
15130de9b25SBarry Smith 
152fb455bf4SPatrick Sanan     Collective
15330de9b25SBarry Smith 
15430de9b25SBarry Smith     Input Parameters:
15530de9b25SBarry Smith +   comm - communicator that shares GUI
15676280437SVaclav Hapla .   prefix - options prefix for all options displayed on window (optional)
15730de9b25SBarry Smith .   title - short descriptive text, for example "Krylov Solver Options"
15887497f52SBarry Smith -   mansec - section of manual pages for options, for example `KSP` (optional)
15930de9b25SBarry Smith 
16030de9b25SBarry Smith     Level: intermediate
16130de9b25SBarry Smith 
16295452b02SPatrick Sanan     Notes:
163d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
164d0609cedSBarry Smith 
16587497f52SBarry Smith     The set of queries needs to be ended by a call to `PetscOptionsEnd()`.
166fb455bf4SPatrick Sanan 
16787497f52SBarry Smith     One can add subheadings with `PetscOptionsHeadBegin()`.
16830de9b25SBarry Smith 
16995452b02SPatrick Sanan     Developer Notes:
17016a05f60SBarry Smith     `PetscOptionsPublish` is set in `PetscOptionsCheckInitial_Private()` with `-saws_options`. When `PetscOptionsPublish` is set the
17116a05f60SBarry Smith     loop between `PetscOptionsBegin()` and `PetscOptionsEnd()` is run THREE times with `PetscOptionsPublishCount` of values -1,0,1.
17216a05f60SBarry Smith      Otherwise the loop is run ONCE with a `PetscOptionsPublishCount` of 1.
17316a05f60SBarry Smith +      \-1 - `PetscOptionsInt()` etc. just call `PetscOptionsGetInt()` etc.
17416a05f60SBarry Smith .      0   - The GUI objects are created in `PetscOptionsInt()` etc. and displayed in `PetscOptionsEnd()` and the options
17516a05f60SBarry Smith               database updated with user changes; `PetscOptionsGetInt()` etc. are also called.
17616a05f60SBarry Smith -      1   - `PetscOptionsInt()` etc. again call `PetscOptionsGetInt()` etc. (possibly getting new values), in addition the help message and
177fb455bf4SPatrick Sanan               default values are printed if -help was given.
17816a05f60SBarry Smith      When `PetscOptionsObject.changedmethod` is set this causes `PetscOptionsPublishCount` to be reset to -2 (so in the next loop iteration it is -1)
17916a05f60SBarry Smith      and the whole process is repeated. This is to handle when, for example, the `KSPType` is changed thus changing the list of
18016a05f60SBarry Smith      options available so they need to be redisplayed so the user can change the. Changing `PetscOptionsObjects.changedmethod` is never
181fb455bf4SPatrick Sanan      currently set.
182aee2cecaSBarry Smith 
1834bb2516aSBarry Smith      Fortran Note:
184feaf08eaSBarry Smith      Returns ierr error code as the final argument per PETSc Fortran API
1854bb2516aSBarry Smith 
186db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
187db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
188db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
189db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
190c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
191db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
192db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()`
19330de9b25SBarry Smith M*/
1949371c9d4SSatish Balay   #define PetscOptionsBegin(comm, prefix, mess, sec) \
1959371c9d4SSatish Balay     do { \
196ce78bad3SBarry Smith       struct _n_PetscOptionItems PetscOptionsObjectBase; \
197ce78bad3SBarry Smith       PetscOptionItems           PetscOptionsObject = &PetscOptionsObjectBase; \
198d0609cedSBarry Smith       PetscCall(PetscMemzero(PetscOptionsObject, sizeof(*PetscOptionsObject))); \
199e55864a3SBarry Smith       for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \
2009566063dSJacob Faibussowitsch         PetscCall(PetscOptionsBegin_Private(PetscOptionsObject, comm, prefix, mess, sec))
20130de9b25SBarry Smith 
2025fefd1ebSJed Brown   /*MC
2035fefd1ebSJed Brown     PetscObjectOptionsBegin - Begins a set of queries on the options database that are related and should be
2045fefd1ebSJed Brown     displayed on the same window of a GUI that allows the user to set the options interactively.
2055fefd1ebSJed Brown 
206f2ba6396SBarry Smith     Synopsis:
207aaa7dc30SBarry Smith     #include <petscoptions.h>
208f2ba6396SBarry Smith     PetscErrorCode PetscObjectOptionsBegin(PetscObject obj)
2095fefd1ebSJed Brown 
210c3339decSBarry Smith     Collective
2115fefd1ebSJed Brown 
2122fe279fdSBarry Smith     Input Parameter:
2135fefd1ebSJed Brown .   obj - object to set options for
2145fefd1ebSJed Brown 
2155fefd1ebSJed Brown     Level: intermediate
2165fefd1ebSJed Brown 
21795452b02SPatrick Sanan     Notes:
218d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
219d0609cedSBarry Smith 
22087497f52SBarry Smith     Needs to be ended by a call the `PetscOptionsEnd()`
221d0609cedSBarry Smith 
22287497f52SBarry Smith     Can add subheadings with `PetscOptionsHeadBegin()`
2235fefd1ebSJed Brown 
224db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
225db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
226db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
227db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
228c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
229db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
230db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2315fefd1ebSJed Brown M*/
2329371c9d4SSatish Balay   #define PetscObjectOptionsBegin(obj) \
2339371c9d4SSatish Balay     do { \
234ce78bad3SBarry Smith       struct _n_PetscOptionItems PetscOptionsObjectBase; \
235ce78bad3SBarry Smith       PetscOptionItems           PetscOptionsObject = &PetscOptionsObjectBase; \
236ce78bad3SBarry Smith       PetscOptionsObject->options                   = ((PetscObject)obj)->options; \
237e55864a3SBarry Smith       for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \
238dbbe0bcdSBarry Smith         PetscCall(PetscObjectOptionsBegin_Private(obj, PetscOptionsObject))
2393194b578SJed Brown 
24030de9b25SBarry Smith   /*MC
24130de9b25SBarry Smith     PetscOptionsEnd - Ends a set of queries on the options database that are related and should be
24230de9b25SBarry Smith      displayed on the same window of a GUI that allows the user to set the options interactively.
24330de9b25SBarry Smith 
244f2ba6396SBarry Smith     Synopsis:
245aaa7dc30SBarry Smith      #include <petscoptions.h>
246f2ba6396SBarry Smith      PetscErrorCode PetscOptionsEnd(void)
24730de9b25SBarry Smith 
2487cdbe19fSJose E. Roman     Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()`
2497cdbe19fSJose E. Roman 
25030de9b25SBarry Smith     Level: intermediate
25130de9b25SBarry Smith 
25295452b02SPatrick Sanan     Notes:
25387497f52SBarry Smith     Needs to be preceded by a call to `PetscOptionsBegin()` or `PetscObjectOptionsBegin()`
25430de9b25SBarry Smith 
255d0609cedSBarry Smith     This is a macro that handles its own error checking, it does not return an error code.
256d0609cedSBarry Smith 
2574bb2516aSBarry Smith     Fortran Note:
258feaf08eaSBarry Smith     Returns ierr error code as the final argument per PETSc Fortran API
2594bb2516aSBarry Smith 
260db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
261db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
262db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
263db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsHeadBegin()`,
264c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
265db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
266db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()`
26730de9b25SBarry Smith M*/
2689371c9d4SSatish Balay   #define PetscOptionsEnd() \
2699371c9d4SSatish Balay     PetscCall(PetscOptionsEnd_Private(PetscOptionsObject)); \
2709371c9d4SSatish Balay     } \
2719371c9d4SSatish Balay     } \
2729371c9d4SSatish Balay     while (0)
2735f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
27430de9b25SBarry Smith 
275ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems, MPI_Comm, const char[], const char[], const char[]);
276ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject, PetscOptionItems);
277ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems);
278ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHeadBegin(PetscOptionItems, const char[]);
27930de9b25SBarry Smith 
2805f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
2819371c9d4SSatish Balay template <typename... T>
2829371c9d4SSatish Balay void PetscOptionsHeadBegin(T...);
283d0609cedSBarry Smith void PetscOptionsHeadEnd(void);
2849371c9d4SSatish Balay template <typename... T>
2859371c9d4SSatish Balay PetscErrorCode PetscOptionsEnum(T...);
2869371c9d4SSatish Balay template <typename... T>
2879371c9d4SSatish Balay PetscErrorCode PetscOptionsInt(T...);
2889371c9d4SSatish Balay template <typename... T>
2899371c9d4SSatish Balay PetscErrorCode PetscOptionsBoundedInt(T...);
2909371c9d4SSatish Balay template <typename... T>
2919371c9d4SSatish Balay PetscErrorCode PetscOptionsRangeInt(T...);
2929371c9d4SSatish Balay template <typename... T>
2939371c9d4SSatish Balay PetscErrorCode PetscOptionsReal(T...);
2949371c9d4SSatish Balay template <typename... T>
2959371c9d4SSatish Balay PetscErrorCode PetscOptionsScalar(T...);
2969371c9d4SSatish Balay template <typename... T>
2979371c9d4SSatish Balay PetscErrorCode PetscOptionsName(T...);
2989371c9d4SSatish Balay template <typename... T>
2999371c9d4SSatish Balay PetscErrorCode PetscOptionsString(T...);
3009371c9d4SSatish Balay template <typename... T>
3019371c9d4SSatish Balay PetscErrorCode PetscOptionsBool(T...);
3029371c9d4SSatish Balay template <typename... T>
3039371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(T...);
3049371c9d4SSatish Balay template <typename... T>
3059371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroup(T...);
3069371c9d4SSatish Balay template <typename... T>
3079371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(T...);
3089371c9d4SSatish Balay template <typename... T>
3099371c9d4SSatish Balay PetscErrorCode PetscOptionsFList(T...);
3109371c9d4SSatish Balay template <typename... T>
3119371c9d4SSatish Balay PetscErrorCode PetscOptionsEList(T...);
3129371c9d4SSatish Balay template <typename... T>
3139371c9d4SSatish Balay PetscErrorCode PetscOptionsRealArray(T...);
3149371c9d4SSatish Balay template <typename... T>
3159371c9d4SSatish Balay PetscErrorCode PetscOptionsScalarArray(T...);
3169371c9d4SSatish Balay template <typename... T>
3179371c9d4SSatish Balay PetscErrorCode PetscOptionsIntArray(T...);
3189371c9d4SSatish Balay template <typename... T>
3199371c9d4SSatish Balay PetscErrorCode PetscOptionsStringArray(T...);
3209371c9d4SSatish Balay template <typename... T>
3219371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolArray(T...);
3229371c9d4SSatish Balay template <typename... T>
3239371c9d4SSatish Balay PetscErrorCode PetscOptionsEnumArray(T...);
3249371c9d4SSatish Balay template <typename... T>
3259371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecated(T...);
3269371c9d4SSatish Balay template <typename... T>
3279371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecatedNoObject(T...);
3285f80ce2aSJacob Faibussowitsch #else
32930de9b25SBarry Smith   /*MC
330d0609cedSBarry Smith    PetscOptionsHeadBegin - Puts a heading before listing any more published options. Used, for example,
33116a05f60SBarry Smith    in `KSPSetFromOptions_GMRES()`.
332d0609cedSBarry Smith 
33316a05f60SBarry Smith    Logically Collective on the communicator passed in `PetscOptionsBegin()`
334d0609cedSBarry Smith 
335d0609cedSBarry Smith    Input Parameter:
336d0609cedSBarry Smith .  head - the heading text
337d0609cedSBarry Smith 
33887497f52SBarry Smith    Level: developer
339d0609cedSBarry Smith 
340d0609cedSBarry Smith    Notes:
341d0609cedSBarry Smith    Handles errors directly, hence does not return an error code
342d0609cedSBarry Smith 
34316a05f60SBarry Smith    Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`, and `PetscOptionsObject` created in `PetscOptionsBegin()` should be the first argument
344d0609cedSBarry Smith 
34587497f52SBarry Smith    Must be followed by a call to `PetscOptionsHeadEnd()` in the same function.
346d0609cedSBarry Smith 
347db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
348db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
349db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
350c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
351db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
352db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
353a54bb2a9SBarry Smith M*/
3549371c9d4SSatish Balay   #define PetscOptionsHeadBegin(PetscOptionsObject, head) \
3559371c9d4SSatish Balay     do { \
35648a46eb9SPierre Jolivet       if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, "  %s\n", head)); \
357d0609cedSBarry Smith     } while (0)
358d0609cedSBarry Smith 
359edd03b47SJacob Faibussowitsch   #define PetscOptionsHead(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadBegin()", ) PetscOptionsHeadBegin(__VA_ARGS__)
360d0609cedSBarry Smith 
361d0609cedSBarry Smith   /*MC
36287497f52SBarry Smith      PetscOptionsHeadEnd - Ends a section of options begun with `PetscOptionsHeadBegin()`
36316a05f60SBarry Smith      See, for example, `KSPSetFromOptions_GMRES()`.
36430de9b25SBarry Smith 
365f2ba6396SBarry Smith      Synopsis:
366aaa7dc30SBarry Smith      #include <petscoptions.h>
367d0609cedSBarry Smith      PetscErrorCode PetscOptionsHeadEnd(void)
36830de9b25SBarry Smith 
3697cdbe19fSJose E. Roman      Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()`
3707cdbe19fSJose E. Roman 
37130de9b25SBarry Smith      Level: intermediate
37230de9b25SBarry Smith 
37395452b02SPatrick Sanan      Notes:
37487497f52SBarry Smith      Must be between a `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` and a `PetscOptionsEnd()`
37530de9b25SBarry Smith 
37687497f52SBarry Smith      Must be preceded by a call to `PetscOptionsHeadBegin()` in the same function.
37730de9b25SBarry Smith 
37887497f52SBarry Smith      This needs to be used only if the code below `PetscOptionsHeadEnd()` can be run ONLY once.
37916a05f60SBarry Smith      See, for example, `PCSetFromOptions_Composite()`. This is a `return(0)` in it for early exit
380b52f573bSBarry Smith      from the function.
381b52f573bSBarry Smith 
38256752e42SBarry Smith      This is only for use with the PETSc options GUI
383b52f573bSBarry Smith 
384db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
385db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
386db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
387c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
388db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
389db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()`
39030de9b25SBarry Smith M*/
3919371c9d4SSatish Balay   #define PetscOptionsHeadEnd() \
3929371c9d4SSatish Balay     do { \
3933ba16761SJacob Faibussowitsch       if (PetscOptionsObject->count != 1) PetscFunctionReturn(PETSC_SUCCESS); \
3949371c9d4SSatish Balay     } while (0)
395d0609cedSBarry Smith 
396edd03b47SJacob Faibussowitsch   #define PetscOptionsTail(...)                                                     PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadEnd()", ) PetscOptionsHeadEnd(__VA_ARGS__)
397186905e3SBarry Smith 
3985305534fSZach Atkins /*MC
3995305534fSZach Atkins   PetscOptionsEnum - Gets the enum value for a particular option in the database.
4005305534fSZach Atkins 
4015305534fSZach Atkins   Synopsis:
4025305534fSZach Atkins   #include <petscoptions.h>
4035305534fSZach Atkins   PetscErrorCode PetscOptionsEnum(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum currentvalue, PetscEnum *value, PetscBool *set)
4045305534fSZach Atkins 
4055305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4065305534fSZach Atkins 
4075305534fSZach Atkins   Input Parameters:
4085305534fSZach Atkins + opt          - option name
4095305534fSZach Atkins . text         - short string that describes the option
4105305534fSZach Atkins . man          - manual page with additional information on option
4115305534fSZach Atkins . list         - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
4125305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
4135305534fSZach Atkins .vb
4145305534fSZach Atkins                  PetscOptionsEnum(..., obj->value,&object->value,...) or
4155305534fSZach Atkins                  value = defaultvalue
4169314d9b7SBarry Smith                  PetscOptionsEnum(..., value,&value,&set);
4179314d9b7SBarry Smith                  if (set) {
4185305534fSZach Atkins .ve
4195305534fSZach Atkins 
4205305534fSZach Atkins   Output Parameters:
4215305534fSZach Atkins + value - the  value to return
4225305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
4235305534fSZach Atkins 
4245305534fSZach Atkins   Level: beginner
4255305534fSZach Atkins 
4265305534fSZach Atkins   Notes:
4275305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
4285305534fSZach Atkins 
4299314d9b7SBarry Smith   `list` is usually something like `PCASMTypes` or some other predefined list of enum names
4305305534fSZach Atkins 
4315305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
4329314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
4335305534fSZach Atkins 
4345305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
4355305534fSZach Atkins 
4365305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
4375305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
4385305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
4395305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
4405305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
4415305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
4425305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
4435305534fSZach Atkins M*/
4441ff8fb82SZach Atkins   #define PetscOptionsEnum(opt, text, man, list, currentvalue, value, set)          PetscOptionsEnum_Private(PetscOptionsObject, opt, text, man, list, currentvalue, value, set)
4455305534fSZach Atkins 
4465305534fSZach Atkins /*MC
4475305534fSZach Atkins   PetscOptionsInt - Gets the integer value for a particular option in the database.
4485305534fSZach Atkins 
4495305534fSZach Atkins   Synopsis:
4505305534fSZach Atkins   #include <petscoptions.h>
4515305534fSZach Atkins   PetscErrorCode PetscOptionsInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set)
4525305534fSZach Atkins 
4535305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4545305534fSZach Atkins 
4555305534fSZach Atkins   Input Parameters:
4565305534fSZach Atkins + opt          - option name
4575305534fSZach Atkins . text         - short string that describes the option
4585305534fSZach Atkins . man          - manual page with additional information on option
4595305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
4605305534fSZach Atkins .vb
4615305534fSZach Atkins                  PetscOptionsInt(..., obj->value, &obj->value, ...) or
4625305534fSZach Atkins                  value = defaultvalue
4639314d9b7SBarry Smith                  PetscOptionsInt(..., value, &value, &set);
4649314d9b7SBarry Smith                  if (set) {
4655305534fSZach Atkins .ve
4665305534fSZach Atkins 
4675305534fSZach Atkins   Output Parameters:
4685305534fSZach Atkins + value - the integer value to return
4695305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
4705305534fSZach Atkins 
4715305534fSZach Atkins   Level: beginner
4725305534fSZach Atkins 
4735305534fSZach Atkins   Notes:
4745305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
4759314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
4765305534fSZach Atkins 
4775305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
4785305534fSZach Atkins 
4795305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
4805305534fSZach Atkins 
4815305534fSZach Atkins .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
4825305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
4835305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
4845305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
4855305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
4865305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
48752ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
4885305534fSZach Atkins M*/
4891690c2aeSBarry Smith   #define PetscOptionsInt(opt, text, man, currentvalue, value, set)                 PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_INT_MIN, PETSC_INT_MAX)
4905305534fSZach Atkins 
4915305534fSZach Atkins /*MC
4926497c311SBarry Smith   PetscOptionsMPIInt - Gets the MPI integer value for a particular option in the database.
4936497c311SBarry Smith 
4946497c311SBarry Smith   Synopsis:
4956497c311SBarry Smith   #include <petscoptions.h>
4966497c311SBarry Smith   PetscErrorCode PetscOptionsMPIInt(const char opt[], const char text[], const char man[], PetscMPIInt currentvalue, PetscMPIInt *value, PetscBool *set)
4976497c311SBarry Smith 
4986497c311SBarry Smith   Logically Collective on the communicator passed in `PetscOptionsBegin()`
4996497c311SBarry Smith 
5006497c311SBarry Smith   Input Parameters:
5016497c311SBarry Smith + opt          - option name
5026497c311SBarry Smith . text         - short string that describes the option
5036497c311SBarry Smith . man          - manual page with additional information on option
5046497c311SBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5056497c311SBarry Smith .vb
5066497c311SBarry Smith                  PetscOptionsInt(..., obj->value, &obj->value, ...) or
5076497c311SBarry Smith                  value = defaultvalue
5086497c311SBarry Smith                  PetscOptionsInt(..., value, &value, &set);
5096497c311SBarry Smith                  if (set) {
5106497c311SBarry Smith .ve
5116497c311SBarry Smith 
5126497c311SBarry Smith   Output Parameters:
5136497c311SBarry Smith + value - the MPI integer value to return
5146497c311SBarry Smith - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
5156497c311SBarry Smith 
5166497c311SBarry Smith   Level: beginner
5176497c311SBarry Smith 
5186497c311SBarry Smith   Notes:
5196497c311SBarry Smith   If the user does not supply the option at all `value` is NOT changed. Thus
5206497c311SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
5216497c311SBarry Smith 
5226497c311SBarry Smith   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
5236497c311SBarry Smith 
5246497c311SBarry Smith   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
5256497c311SBarry Smith 
5266497c311SBarry Smith .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
5276497c311SBarry Smith           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
5286497c311SBarry Smith           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
5296497c311SBarry Smith           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
5306497c311SBarry Smith           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
5316497c311SBarry Smith           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
5326497c311SBarry Smith           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
5336497c311SBarry Smith M*/
5346497c311SBarry 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)
5356497c311SBarry Smith 
5366497c311SBarry Smith /*MC
53752ce0ab5SPierre Jolivet    PetscOptionsBoundedInt - Gets an integer value greater than or equal to a given bound for a particular option in the database.
5385305534fSZach Atkins 
5395305534fSZach Atkins    Synopsis:
5405305534fSZach Atkins    #include <petscoptions.h>
5419314d9b7SBarry Smith    PetscErrorCode  PetscOptionsBoundedInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt bound)
5425305534fSZach Atkins 
5435305534fSZach Atkins    Logically Collective on the communicator passed in `PetscOptionsBegin()`
5445305534fSZach Atkins 
5455305534fSZach Atkins    Input Parameters:
5465305534fSZach Atkins +  opt          - option name
5475305534fSZach Atkins .  text         - short string that describes the option
5485305534fSZach Atkins .  man          - manual page with additional information on option
5495305534fSZach Atkins .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5505305534fSZach Atkins .vb
55152ce0ab5SPierre Jolivet   PetscOptionsBoundedInt(..., obj->value, &obj->value, ...)
5525305534fSZach Atkins .ve
5535305534fSZach Atkins or
5545305534fSZach Atkins .vb
5555305534fSZach Atkins   value = defaultvalue
55652ce0ab5SPierre Jolivet   PetscOptionsBoundedInt(..., value, &value, &set, ...);
5579314d9b7SBarry Smith   if (set) {
5585305534fSZach Atkins .ve
55952ce0ab5SPierre Jolivet -  bound - the requested value should be greater than or equal to this bound or an error is generated
5605305534fSZach Atkins 
5615305534fSZach Atkins    Output Parameters:
5625305534fSZach Atkins +  value - the integer value to return
5639314d9b7SBarry Smith -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
5645305534fSZach Atkins 
5655305534fSZach Atkins    Level: beginner
5665305534fSZach Atkins 
5675305534fSZach Atkins    Notes:
5685305534fSZach Atkins    If the user does not supply the option at all `value` is NOT changed. Thus
5699314d9b7SBarry Smith    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
5705305534fSZach Atkins 
5715305534fSZach Atkins    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
5725305534fSZach Atkins 
573af27ebaaSBarry Smith    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
5745305534fSZach Atkins 
5755305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
5765305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
5775305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
5785305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
5795305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
5805305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
58152ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
5825305534fSZach Atkins M*/
5831690c2aeSBarry Smith   #define PetscOptionsBoundedInt(opt, text, man, currentvalue, value, set, lb)      PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_INT_MAX)
5845305534fSZach Atkins 
5855305534fSZach Atkins /*MC
5865305534fSZach Atkins    PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database.
5875305534fSZach Atkins 
5885305534fSZach Atkins    Synopsis:
5895305534fSZach Atkins    #include <petscoptions.h>
5909314d9b7SBarry Smith    PetscErrorCode PetscOptionsRangeInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt lb, PetscInt ub)
5915305534fSZach Atkins 
5925305534fSZach Atkins    Logically Collective on the communicator passed in `PetscOptionsBegin()`
5935305534fSZach Atkins 
5945305534fSZach Atkins    Input Parameters:
5955305534fSZach Atkins +  opt          - option name
5965305534fSZach Atkins .  text         - short string that describes the option
5975305534fSZach Atkins .  man          - manual page with additional information on option
5985305534fSZach Atkins .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
5995305534fSZach Atkins .vb
60052ce0ab5SPierre Jolivet   PetscOptionsRangeInt(..., obj->value, &obj->value, ...)
60152ce0ab5SPierre Jolivet .ve
60252ce0ab5SPierre Jolivet or
60352ce0ab5SPierre Jolivet .vb
6045305534fSZach Atkins   value = defaultvalue
60552ce0ab5SPierre Jolivet   PetscOptionsRangeInt(..., value, &value, &set, ...);
6069314d9b7SBarry Smith   if (set) {
6075305534fSZach Atkins .ve
6085305534fSZach Atkins .  lb - the lower bound, provided value must be greater than or equal to this value or an error is generated
6095305534fSZach Atkins -  ub - the upper bound, provided value must be less than or equal to this value or an error is generated
6105305534fSZach Atkins 
6115305534fSZach Atkins    Output Parameters:
6125305534fSZach Atkins +  value - the integer value to return
6139314d9b7SBarry Smith -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
6145305534fSZach Atkins 
6155305534fSZach Atkins    Level: beginner
6165305534fSZach Atkins 
6175305534fSZach Atkins    Notes:
6185305534fSZach Atkins    If the user does not supply the option at all `value` is NOT changed. Thus
6199314d9b7SBarry Smith    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
6205305534fSZach Atkins 
6215305534fSZach Atkins    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
6225305534fSZach Atkins 
623af27ebaaSBarry Smith    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
6245305534fSZach Atkins 
6255305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
6265305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()`
6275305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
6285305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
6295305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
6305305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
63152ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
6325305534fSZach Atkins M*/
6331ff8fb82SZach Atkins   #define PetscOptionsRangeInt(opt, text, man, currentvalue, value, set, lb, ub)    PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub)
6345305534fSZach Atkins 
6355305534fSZach Atkins /*MC
63652ce0ab5SPierre Jolivet   PetscOptionsReal - Gets a `PetscReal` value for a particular option in the database.
6375305534fSZach Atkins 
6385305534fSZach Atkins   Synopsis:
6395305534fSZach Atkins   #include <petscoptions.h>
6405305534fSZach Atkins   PetscErrorCode PetscOptionsReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set)
6415305534fSZach Atkins 
6425305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
6435305534fSZach Atkins 
6445305534fSZach Atkins   Input Parameters:
6455305534fSZach Atkins + opt          - option name
6465305534fSZach Atkins . text         - short string that describes the option
6475305534fSZach Atkins . man          - manual page with additional information on option
6485305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6495305534fSZach Atkins .vb
6505305534fSZach Atkins                  PetscOptionsReal(..., obj->value,&obj->value,...) or
6515305534fSZach Atkins                  value = defaultvalue
6529314d9b7SBarry Smith                  PetscOptionsReal(..., value,&value,&set);
6539314d9b7SBarry Smith                  if (set) {
6545305534fSZach Atkins .ve
6555305534fSZach Atkins 
6565305534fSZach Atkins   Output Parameters:
6575305534fSZach Atkins + value - the value to return
6585305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
6595305534fSZach Atkins 
6605305534fSZach Atkins   Level: beginner
6615305534fSZach Atkins 
6625305534fSZach Atkins   Notes:
6635305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
6649314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
6655305534fSZach Atkins 
6665305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
6675305534fSZach Atkins 
668af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
6695305534fSZach Atkins 
6705305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
6715305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
6725305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
6735305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
6745305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
6755305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
67652ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()`
6775305534fSZach Atkins M*/
67852ce0ab5SPierre Jolivet   #define PetscOptionsReal(opt, text, man, currentvalue, value, set)                PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MIN_REAL, PETSC_MAX_REAL)
67952ce0ab5SPierre Jolivet 
68052ce0ab5SPierre Jolivet /*MC
68152ce0ab5SPierre Jolivet    PetscOptionsBoundedReal - Gets a `PetscReal` value greater than or equal to a given bound for a particular option in the database.
68252ce0ab5SPierre Jolivet 
68352ce0ab5SPierre Jolivet    Synopsis:
68452ce0ab5SPierre Jolivet    #include <petscoptions.h>
68552ce0ab5SPierre Jolivet    PetscErrorCode  PetscOptionsBoundedReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal bound)
68652ce0ab5SPierre Jolivet 
68752ce0ab5SPierre Jolivet    Logically Collective on the communicator passed in `PetscOptionsBegin()`
68852ce0ab5SPierre Jolivet 
68952ce0ab5SPierre Jolivet    Input Parameters:
69052ce0ab5SPierre Jolivet +  opt          - option name
69152ce0ab5SPierre Jolivet .  text         - short string that describes the option
69252ce0ab5SPierre Jolivet .  man          - manual page with additional information on option
69352ce0ab5SPierre Jolivet .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
69452ce0ab5SPierre Jolivet .vb
69552ce0ab5SPierre Jolivet   PetscOptionsBoundedReal(..., obj->value, &obj->value, ...)
69652ce0ab5SPierre Jolivet .ve
69752ce0ab5SPierre Jolivet or
69852ce0ab5SPierre Jolivet .vb
69952ce0ab5SPierre Jolivet   value = defaultvalue
70052ce0ab5SPierre Jolivet   PetscOptionsBoundedReal(..., value, &value, &set, ...);
70152ce0ab5SPierre Jolivet   if (set) {
70252ce0ab5SPierre Jolivet .ve
70352ce0ab5SPierre Jolivet -  bound - the requested value should be greater than or equal to this bound or an error is generated
70452ce0ab5SPierre Jolivet 
70552ce0ab5SPierre Jolivet    Output Parameters:
70652ce0ab5SPierre Jolivet +  value - the real value to return
70752ce0ab5SPierre Jolivet -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
70852ce0ab5SPierre Jolivet 
70952ce0ab5SPierre Jolivet    Level: beginner
71052ce0ab5SPierre Jolivet 
71152ce0ab5SPierre Jolivet    Notes:
71252ce0ab5SPierre Jolivet    If the user does not supply the option at all `value` is NOT changed. Thus
71352ce0ab5SPierre Jolivet    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
71452ce0ab5SPierre Jolivet 
71552ce0ab5SPierre Jolivet    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
71652ce0ab5SPierre Jolivet 
71752ce0ab5SPierre Jolivet    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
71852ce0ab5SPierre Jolivet 
71952ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
72052ce0ab5SPierre Jolivet           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()`
72152ce0ab5SPierre Jolivet           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
72252ce0ab5SPierre Jolivet           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
72352ce0ab5SPierre Jolivet           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
72452ce0ab5SPierre Jolivet           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
72552ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedInt()`, `PetscOptionsRangeReal()`
72652ce0ab5SPierre Jolivet M*/
72752ce0ab5SPierre Jolivet   #define PetscOptionsBoundedReal(opt, text, man, currentvalue, value, set, lb)     PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_MAX_REAL)
72852ce0ab5SPierre Jolivet 
72952ce0ab5SPierre Jolivet /*MC
73052ce0ab5SPierre Jolivet    PetscOptionsRangeReal - Gets a `PetscReal` value within a range of values for a particular option in the database.
73152ce0ab5SPierre Jolivet 
73252ce0ab5SPierre Jolivet    Synopsis:
73352ce0ab5SPierre Jolivet    #include <petscoptions.h>
73452ce0ab5SPierre Jolivet    PetscErrorCode PetscOptionsRangeReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal lb, PetscReal ub)
73552ce0ab5SPierre Jolivet 
73652ce0ab5SPierre Jolivet    Logically Collective on the communicator passed in `PetscOptionsBegin()`
73752ce0ab5SPierre Jolivet 
73852ce0ab5SPierre Jolivet    Input Parameters:
73952ce0ab5SPierre Jolivet +  opt          - option name
74052ce0ab5SPierre Jolivet .  text         - short string that describes the option
74152ce0ab5SPierre Jolivet .  man          - manual page with additional information on option
74252ce0ab5SPierre Jolivet .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
74352ce0ab5SPierre Jolivet .vb
74452ce0ab5SPierre Jolivet   PetscOptionsRangeReal(..., obj->value, &obj->value, ...)
74552ce0ab5SPierre Jolivet .ve
74652ce0ab5SPierre Jolivet or
74752ce0ab5SPierre Jolivet .vb
74852ce0ab5SPierre Jolivet   value = defaultvalue
74952ce0ab5SPierre Jolivet   PetscOptionsRangeReal(..., value, &value, &set, ...);
75052ce0ab5SPierre Jolivet   if (set) {
75152ce0ab5SPierre Jolivet .ve
75252ce0ab5SPierre Jolivet .  lb - the lower bound, provided value must be greater than or equal to this value or an error is generated
75352ce0ab5SPierre Jolivet -  ub - the upper bound, provided value must be less than or equal to this value or an error is generated
75452ce0ab5SPierre Jolivet 
75552ce0ab5SPierre Jolivet    Output Parameters:
75652ce0ab5SPierre Jolivet +  value - the value to return
75752ce0ab5SPierre Jolivet -  set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
75852ce0ab5SPierre Jolivet 
75952ce0ab5SPierre Jolivet    Level: beginner
76052ce0ab5SPierre Jolivet 
76152ce0ab5SPierre Jolivet    Notes:
76252ce0ab5SPierre Jolivet    If the user does not supply the option at all `value` is NOT changed. Thus
76352ce0ab5SPierre Jolivet    you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
76452ce0ab5SPierre Jolivet 
76552ce0ab5SPierre Jolivet    The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
76652ce0ab5SPierre Jolivet 
76752ce0ab5SPierre Jolivet    Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
76852ce0ab5SPierre Jolivet 
76952ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
77052ce0ab5SPierre Jolivet           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()`
77152ce0ab5SPierre Jolivet           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
77252ce0ab5SPierre Jolivet           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
77352ce0ab5SPierre Jolivet           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
77452ce0ab5SPierre Jolivet           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
77552ce0ab5SPierre Jolivet           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRangeInt()`, `PetscOptionsBoundedReal()`
77652ce0ab5SPierre Jolivet M*/
77752ce0ab5SPierre Jolivet   #define PetscOptionsRangeReal(opt, text, man, currentvalue, value, set, lb, ub)   PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub)
7785305534fSZach Atkins 
7795305534fSZach Atkins /*MC
7805305534fSZach Atkins   PetscOptionsScalar - Gets the `PetscScalar` value for a particular option in the database.
7815305534fSZach Atkins 
7825305534fSZach Atkins   Synopsis:
7835305534fSZach Atkins   #include <petscoptions.h>
7845305534fSZach Atkins   PetscErrorCode PetscOptionsScalar(const char opt[], const char text[], const char man[], PetscScalar currentvalue, PetscScalar *value, PetscBool *set)
7855305534fSZach Atkins 
7865305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
7875305534fSZach Atkins 
7885305534fSZach Atkins   Input Parameters:
7895305534fSZach Atkins + opt          - option name
7905305534fSZach Atkins . text         - short string that describes the option
7915305534fSZach Atkins . man          - manual page with additional information on option
7925305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7935305534fSZach Atkins .vb
7945305534fSZach Atkins                  PetscOptionsScalar(..., obj->value,&obj->value,...) or
7955305534fSZach Atkins                  value = defaultvalue
7969314d9b7SBarry Smith                  PetscOptionsScalar(..., value,&value,&set);
7979314d9b7SBarry Smith                  if (set) {
7985305534fSZach Atkins .ve
7995305534fSZach Atkins 
8005305534fSZach Atkins   Output Parameters:
8015305534fSZach Atkins + value - the value to return
8025305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
8035305534fSZach Atkins 
8045305534fSZach Atkins   Level: beginner
8055305534fSZach Atkins 
8065305534fSZach Atkins   Notes:
8075305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
8089314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
8095305534fSZach Atkins 
8105305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
8115305534fSZach Atkins 
812af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8135305534fSZach Atkins 
8145305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8155305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8165305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8175305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8185305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8195305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8205305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8215305534fSZach Atkins M*/
8221ff8fb82SZach Atkins   #define PetscOptionsScalar(opt, text, man, currentvalue, value, set)              PetscOptionsScalar_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
8235305534fSZach Atkins 
8245305534fSZach Atkins /*MC
8255305534fSZach 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
8265305534fSZach Atkins   its value is set to false.
8275305534fSZach Atkins 
8285305534fSZach Atkins   Synopsis:
8295305534fSZach Atkins   #include <petscoptions.h>
8309314d9b7SBarry Smith   PetscErrorCode PetscOptionsName(const char opt[], const char text[], const char man[], PetscBool *set)
8315305534fSZach Atkins 
8325305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
8335305534fSZach Atkins 
8345305534fSZach Atkins   Input Parameters:
8355305534fSZach Atkins + opt  - option name
8365305534fSZach Atkins . text - short string that describes the option
8375305534fSZach Atkins - man  - manual page with additional information on option
8385305534fSZach Atkins 
8395305534fSZach Atkins   Output Parameter:
8409314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
8415305534fSZach Atkins 
8425305534fSZach Atkins   Level: beginner
8435305534fSZach Atkins 
8445305534fSZach Atkins   Note:
845af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8465305534fSZach Atkins 
8475305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8485305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8495305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8505305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8515305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8525305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8535305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8545305534fSZach Atkins M*/
8559314d9b7SBarry Smith   #define PetscOptionsName(opt, text, man, set)                                     PetscOptionsName_Private(PetscOptionsObject, opt, text, man, set)
8565305534fSZach Atkins 
8575305534fSZach Atkins /*MC
8585305534fSZach Atkins   PetscOptionsString - Gets the string value for a particular option in the database.
8595305534fSZach Atkins 
8605305534fSZach Atkins   Synopsis:
8615305534fSZach Atkins   #include <petscoptions.h>
8625305534fSZach Atkins   PetscErrorCode PetscOptionsString(const char opt[], const char text[], const char man[], const char currentvalue[], char value[], size_t len, PetscBool *set)
8635305534fSZach Atkins 
8645305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
8655305534fSZach Atkins 
8665305534fSZach Atkins   Input Parameters:
8675305534fSZach Atkins + opt          - option name
8685305534fSZach Atkins . text         - short string that describes the option
8695305534fSZach Atkins . man          - manual page with additional information on option
8705305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
8715305534fSZach Atkins - len          - length of the result string including null terminator
8725305534fSZach Atkins 
8735305534fSZach Atkins   Output Parameters:
8745305534fSZach Atkins + value - the value to return
8755305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
8765305534fSZach Atkins 
8775305534fSZach Atkins   Level: beginner
8785305534fSZach Atkins 
8795305534fSZach Atkins   Notes:
880af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
8815305534fSZach Atkins 
8829314d9b7SBarry Smith   If the user provided no string (for example `-optionname` `-someotheroption`) `set` is set to `PETSC_TRUE` (and the string is filled with nulls).
8835305534fSZach Atkins 
8845305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
8859314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`.
8865305534fSZach Atkins 
8875305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
8885305534fSZach Atkins 
8895305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
8905305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
8915305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
8925305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
8935305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
8945305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
8955305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
8965305534fSZach Atkins M*/
8971ff8fb82SZach Atkins   #define PetscOptionsString(opt, text, man, currentvalue, value, len, set)         PetscOptionsString_Private(PetscOptionsObject, opt, text, man, currentvalue, value, len, set)
8985305534fSZach Atkins 
8995305534fSZach Atkins /*MC
9005305534fSZach Atkins   PetscOptionsBool - Determines if a particular option is in the database with a true or false
9015305534fSZach Atkins 
9025305534fSZach Atkins   Synopsis:
9035305534fSZach Atkins   #include <petscoptions.h>
9045305534fSZach Atkins   PetscErrorCode PetscOptionsBool(const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool *flg, PetscBool *set)
9055305534fSZach Atkins 
9065305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
9075305534fSZach Atkins 
9085305534fSZach Atkins   Input Parameters:
9095305534fSZach Atkins + opt          - option name
9105305534fSZach Atkins . text         - short string that describes the option
9115305534fSZach Atkins . man          - manual page with additional information on option
9125305534fSZach Atkins - currentvalue - the current value
9135305534fSZach Atkins 
9145305534fSZach Atkins   Output Parameters:
9155305534fSZach Atkins + flg - `PETSC_TRUE` or `PETSC_FALSE`
916*4d81f786SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`, pass `NULL` if not needed
9175305534fSZach Atkins 
9185305534fSZach Atkins   Level: beginner
9195305534fSZach Atkins 
9205305534fSZach Atkins   Notes:
9215305534fSZach Atkins   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
9225305534fSZach Atkins   FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
9235305534fSZach Atkins 
9249314d9b7SBarry 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`
9255305534fSZach Atkins   is equivalent to `-requested_bool true`
9265305534fSZach Atkins 
9275305534fSZach Atkins   If the user does not supply the option at all `flg` is NOT changed. Thus
9289314d9b7SBarry Smith   you should ALWAYS initialize the `flg` variable if you access it without first checking that the `set` flag is `PETSC_TRUE`.
9295305534fSZach Atkins 
9305305534fSZach Atkins   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
9315305534fSZach Atkins 
9325305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
9335305534fSZach Atkins           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
9345305534fSZach Atkins           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
9355305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
9365305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
9375305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
9385305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
9395305534fSZach Atkins M*/
9401ff8fb82SZach Atkins   #define PetscOptionsBool(opt, text, man, currentvalue, value, set)                PetscOptionsBool_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
9415305534fSZach Atkins 
9425305534fSZach Atkins /*MC
943*4d81f786SBarry Smith   PetscOptionsBool3 - Determines if a particular option is in the database with a true, false, or unknown
944*4d81f786SBarry Smith 
945*4d81f786SBarry Smith   Synopsis:
946*4d81f786SBarry Smith   #include <petscoptions.h>
947*4d81f786SBarry Smith   PetscErrorCode PetscOptionsBool3(const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool3 *flg, PetscBool *set)
948*4d81f786SBarry Smith 
949*4d81f786SBarry Smith   Logically Collective on the communicator passed in `PetscOptionsBegin()`
950*4d81f786SBarry Smith 
951*4d81f786SBarry Smith   Input Parameters:
952*4d81f786SBarry Smith + opt          - option name
953*4d81f786SBarry Smith . text         - short string that describes the option
954*4d81f786SBarry Smith . man          - manual page with additional information on option
955*4d81f786SBarry Smith - currentvalue - the current value
956*4d81f786SBarry Smith 
957*4d81f786SBarry Smith   Output Parameters:
958*4d81f786SBarry Smith + flg - `PETSC_BOOL3_TRUE`, `PETSC_BOOL3_FALSE`, or `PETSC_BOOL3_UNKNOWN`
959*4d81f786SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
960*4d81f786SBarry Smith 
961*4d81f786SBarry Smith   Level: beginner
962*4d81f786SBarry Smith 
963*4d81f786SBarry Smith   Notes:
964*4d81f786SBarry Smith   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
965*4d81f786SBarry Smith   FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
966*4d81f786SBarry Smith 
967*4d81f786SBarry Smith   If the option is given, but no value is provided, then `flg` and `set` are both given the value `PETSC_BOOL3_TRUE`. That is `-requested_bool`
968*4d81f786SBarry Smith   is equivalent to `-requested_bool true`
969*4d81f786SBarry Smith 
970*4d81f786SBarry Smith   If the user does not supply the option at all `flg` is NOT changed. Thus
971*4d81f786SBarry Smith   you should ALWAYS initialize the `flg` variable if you access it without first checking that the `set` flag is `PETSC_TRUE`.
972*4d81f786SBarry Smith 
973*4d81f786SBarry Smith   Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
974*4d81f786SBarry Smith 
975*4d81f786SBarry Smith .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
976*4d81f786SBarry Smith           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
977*4d81f786SBarry Smith           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
978*4d81f786SBarry Smith           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
979*4d81f786SBarry Smith           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
980*4d81f786SBarry Smith           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
981*4d81f786SBarry Smith           `PetscOptionsFList()`, `PetscOptionsEList()`
982*4d81f786SBarry Smith M*/
983*4d81f786SBarry Smith   #define PetscOptionsBool3(opt, text, man, currentvalue, value, set)               PetscOptionsBool3_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set)
984*4d81f786SBarry Smith 
985*4d81f786SBarry Smith /*MC
9865305534fSZach Atkins   PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
9875305534fSZach Atkins   which at most a single value can be true.
9885305534fSZach Atkins 
9895305534fSZach Atkins   Synopsis:
9905305534fSZach Atkins   #include <petscoptions.h>
9919314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[], const char text[], const char man[], PetscBool *set)
9925305534fSZach Atkins 
9935305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
9945305534fSZach Atkins 
9955305534fSZach Atkins   Input Parameters:
9965305534fSZach Atkins + opt  - option name
9975305534fSZach Atkins . text - short string that describes the option
9985305534fSZach Atkins - man  - manual page with additional information on option
9995305534fSZach Atkins 
10005305534fSZach Atkins   Output Parameter:
10019314d9b7SBarry Smith . set - whether that option was set or not
10025305534fSZach Atkins 
10035305534fSZach Atkins   Level: intermediate
10045305534fSZach Atkins 
10055305534fSZach Atkins   Notes:
1006af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
10075305534fSZach Atkins 
10085305534fSZach Atkins   Must be followed by 0 or more `PetscOptionsBoolGroup()`s and `PetscOptionsBoolGroupEnd()`
10095305534fSZach Atkins 
10105305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10115305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10125305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10135305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10145305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10155305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
10165305534fSZach Atkins M*/
10179314d9b7SBarry Smith   #define PetscOptionsBoolGroupBegin(opt, text, man, set)                           PetscOptionsBoolGroupBegin_Private(PetscOptionsObject, opt, text, man, set)
10185305534fSZach Atkins 
10195305534fSZach Atkins /*MC
10205305534fSZach Atkins   PetscOptionsBoolGroup - One in a series of logical queries on the options database for
10215305534fSZach Atkins   which at most a single value can be true.
10225305534fSZach Atkins 
10235305534fSZach Atkins   Synopsis:
10245305534fSZach Atkins   #include <petscoptions.h>
10259314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroup(const char opt[], const char text[], const char man[], PetscBool *set)
10265305534fSZach Atkins 
10275305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
10285305534fSZach Atkins 
10295305534fSZach Atkins   Input Parameters:
10305305534fSZach Atkins + opt  - option name
10315305534fSZach Atkins . text - short string that describes the option
10325305534fSZach Atkins - man  - manual page with additional information on option
10335305534fSZach Atkins 
10345305534fSZach Atkins   Output Parameter:
10359314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
10365305534fSZach Atkins 
10375305534fSZach Atkins   Level: intermediate
10385305534fSZach Atkins 
10395305534fSZach Atkins   Notes:
1040af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
10415305534fSZach Atkins 
10425305534fSZach Atkins   Must follow a `PetscOptionsBoolGroupBegin()` and preceded a `PetscOptionsBoolGroupEnd()`
10435305534fSZach Atkins 
10445305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10455305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10465305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10475305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10485305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10495305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
10505305534fSZach Atkins M*/
10519314d9b7SBarry Smith   #define PetscOptionsBoolGroup(opt, text, man, set)                                PetscOptionsBoolGroup_Private(PetscOptionsObject, opt, text, man, set)
10525305534fSZach Atkins 
10535305534fSZach Atkins /*MC
10545305534fSZach Atkins   PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
10555305534fSZach Atkins   which at most a single value can be true.
10565305534fSZach Atkins 
10575305534fSZach Atkins   Synopsis:
10585305534fSZach Atkins   #include <petscoptions.h>
10599314d9b7SBarry Smith   PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[], const char text[], const char man[], PetscBool  *set)
10605305534fSZach Atkins 
10615305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
10625305534fSZach Atkins 
10635305534fSZach Atkins   Input Parameters:
10645305534fSZach Atkins + opt  - option name
10655305534fSZach Atkins . text - short string that describes the option
10665305534fSZach Atkins - man  - manual page with additional information on option
10675305534fSZach Atkins 
10685305534fSZach Atkins   Output Parameter:
10699314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE`
10705305534fSZach Atkins 
10715305534fSZach Atkins   Level: intermediate
10725305534fSZach Atkins 
10735305534fSZach Atkins   Notes:
1074af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
10755305534fSZach Atkins 
10765305534fSZach Atkins   Must follow a `PetscOptionsBoolGroupBegin()`
10775305534fSZach Atkins 
10785305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
10795305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
10805305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
10815305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
10825305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
10835305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
10845305534fSZach Atkins M*/
10859314d9b7SBarry Smith   #define PetscOptionsBoolGroupEnd(opt, text, man, set)                             PetscOptionsBoolGroupEnd_Private(PetscOptionsObject, opt, text, man, set)
10865305534fSZach Atkins 
10875305534fSZach Atkins /*MC
10885305534fSZach Atkins   PetscOptionsFList - Puts a list of option values that a single one may be selected from
10895305534fSZach Atkins 
10905305534fSZach Atkins   Synopsis:
10915305534fSZach Atkins   #include <petscoptions.h>
10925305534fSZach Atkins   PetscErrorCode PetscOptionsFList(const char opt[], const char ltext[], const char man[], PetscFunctionList list, const char currentvalue[], char value[], size_t len, PetscBool *set)
10935305534fSZach Atkins 
10945305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
10955305534fSZach Atkins 
10965305534fSZach Atkins   Input Parameters:
10975305534fSZach Atkins + opt          - option name
10985305534fSZach Atkins . ltext        - short string that describes the option
10995305534fSZach Atkins . man          - manual page with additional information on option
11005305534fSZach Atkins . list         - the possible choices
11015305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11025305534fSZach Atkins .vb
11039314d9b7SBarry Smith                  PetscOptionsFlist(..., obj->value,value,len,&set);
11049314d9b7SBarry Smith                  if (set) {
11055305534fSZach Atkins .ve
11065305534fSZach Atkins - len          - the length of the character array value
11075305534fSZach Atkins 
11085305534fSZach Atkins   Output Parameters:
11095305534fSZach Atkins + value - the value to return
11105305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
11115305534fSZach Atkins 
11125305534fSZach Atkins   Level: intermediate
11135305534fSZach Atkins 
11145305534fSZach Atkins   Notes:
1115af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
11165305534fSZach Atkins 
11175305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
11189314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`.
11195305534fSZach Atkins 
11205305534fSZach Atkins   The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically.
11215305534fSZach Atkins 
11225305534fSZach Atkins   See `PetscOptionsEList()` for when the choices are given in a string array
11235305534fSZach Atkins 
11245305534fSZach Atkins   To get a listing of all currently specified options,
11255305534fSZach Atkins   see `PetscOptionsView()` or `PetscOptionsGetAll()`
11265305534fSZach Atkins 
1127af27ebaaSBarry Smith   Developer Note:
11285305534fSZach Atkins   This cannot check for invalid selection because of things like `MATAIJ` that are not included in the list
11295305534fSZach Atkins 
11305305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
11315305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
11325305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
11335305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
11345305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
11355305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()`
11365305534fSZach Atkins M*/
11371ff8fb82SZach Atkins   #define PetscOptionsFList(opt, ltext, man, list, currentvalue, value, len, set)   PetscOptionsFList_Private(PetscOptionsObject, opt, ltext, man, list, currentvalue, value, len, set)
11385305534fSZach Atkins 
11395305534fSZach Atkins /*MC
11405305534fSZach Atkins   PetscOptionsEList - Puts a list of option values that a single one may be selected from
11415305534fSZach Atkins 
11425305534fSZach Atkins   Synopsis:
11435305534fSZach Atkins   #include <petscoptions.h>
11445305534fSZach 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)
11455305534fSZach Atkins 
11465305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
11475305534fSZach Atkins 
11485305534fSZach Atkins   Input Parameters:
11495305534fSZach Atkins + opt          - option name
11505305534fSZach Atkins . ltext        - short string that describes the option
11515305534fSZach Atkins . man          - manual page with additional information on option
11525305534fSZach Atkins . list         - the possible choices (one of these must be selected, anything else is invalid)
11535305534fSZach Atkins . ntext        - number of choices
11545305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
11555305534fSZach Atkins .vb
11569314d9b7SBarry Smith                  PetscOptionsEList(..., obj->value,&value,&set);
11579314d9b7SBarry Smith .ve                 if (set) {
11585305534fSZach Atkins 
11595305534fSZach Atkins   Output Parameters:
11605305534fSZach Atkins + value - the index of the value to return
11615305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
11625305534fSZach Atkins 
11635305534fSZach Atkins   Level: intermediate
11645305534fSZach Atkins 
11655305534fSZach Atkins   Notes:
1166af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
11675305534fSZach Atkins 
11685305534fSZach Atkins   If the user does not supply the option at all `value` is NOT changed. Thus
11699314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`.
11705305534fSZach Atkins 
11715305534fSZach Atkins   See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList()`
11725305534fSZach Atkins 
11735305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
11745305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
11755305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
11765305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
11775305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
11785305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEnum()`
11795305534fSZach Atkins M*/
11801ff8fb82SZach Atkins   #define PetscOptionsEList(opt, ltext, man, list, ntext, currentvalue, value, set) PetscOptionsEList_Private(PetscOptionsObject, opt, ltext, man, list, ntext, currentvalue, value, set)
11815305534fSZach Atkins 
11825305534fSZach Atkins /*MC
11835305534fSZach Atkins   PetscOptionsRealArray - Gets an array of double values for a particular
11845305534fSZach Atkins   option in the database. The values must be separated with commas with
11855305534fSZach Atkins   no intervening spaces.
11865305534fSZach Atkins 
11875305534fSZach Atkins   Synopsis:
11885305534fSZach Atkins   #include <petscoptions.h>
11895305534fSZach Atkins   PetscErrorCode PetscOptionsRealArray(const char opt[], const char text[], const char man[], PetscReal value[], PetscInt *n, PetscBool *set)
11905305534fSZach Atkins 
11915305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
11925305534fSZach Atkins 
11935305534fSZach Atkins   Input Parameters:
11945305534fSZach Atkins + opt  - the option one is seeking
11955305534fSZach Atkins . text - short string describing option
11965305534fSZach Atkins . man  - manual page for option
11975305534fSZach Atkins - n    - maximum number of values that value has room for
11985305534fSZach Atkins 
11995305534fSZach Atkins   Output Parameters:
12005305534fSZach Atkins + value - location to copy values
12015305534fSZach Atkins . n     - actual number of values found
12025305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
12035305534fSZach Atkins 
12045305534fSZach Atkins   Level: beginner
12055305534fSZach Atkins 
12065305534fSZach Atkins   Note:
1207af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
12085305534fSZach Atkins 
12095305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12105305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12115305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12125305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12135305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12145305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12155305534fSZach Atkins M*/
1216c857d48bSBarry Smith   #define PetscOptionsRealArray(opt, text, man, value, n, set)                      PetscOptionsRealArray_Private(PetscOptionsObject, opt, text, man, value, n, set)
12175305534fSZach Atkins 
12185305534fSZach Atkins /*MC
12195305534fSZach Atkins   PetscOptionsScalarArray - Gets an array of `PetscScalar` values for a particular
12205305534fSZach Atkins   option in the database. The values must be separated with commas with
12215305534fSZach Atkins   no intervening spaces.
12225305534fSZach Atkins 
12235305534fSZach Atkins   Synopsis:
12245305534fSZach Atkins   #include <petscoptions.h>
12255305534fSZach Atkins   PetscErrorCode PetscOptionsScalarArray(const char opt[], const char text[], const char man[], PetscScalar value[], PetscInt *n, PetscBool *set)
12265305534fSZach Atkins 
12275305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
12285305534fSZach Atkins 
12295305534fSZach Atkins   Input Parameters:
12305305534fSZach Atkins + opt  - the option one is seeking
12315305534fSZach Atkins . text - short string describing option
12325305534fSZach Atkins . man  - manual page for option
12335305534fSZach Atkins - n    - maximum number of values allowed in the value array
12345305534fSZach Atkins 
12355305534fSZach Atkins   Output Parameters:
12365305534fSZach Atkins + value - location to copy values
12375305534fSZach Atkins . n     - actual number of values found
12385305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
12395305534fSZach Atkins 
12405305534fSZach Atkins   Level: beginner
12415305534fSZach Atkins 
12425305534fSZach Atkins   Note:
1243af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
12445305534fSZach Atkins 
12455305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12465305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12475305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12485305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12495305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12505305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12515305534fSZach Atkins M*/
1252c857d48bSBarry Smith   #define PetscOptionsScalarArray(opt, text, man, value, n, set)                    PetscOptionsScalarArray_Private(PetscOptionsObject, opt, text, man, value, n, set)
12535305534fSZach Atkins 
12545305534fSZach Atkins /*MC
12555305534fSZach Atkins   PetscOptionsIntArray - Gets an array of integers for a particular
12565305534fSZach Atkins   option in the database.
12575305534fSZach Atkins 
12585305534fSZach Atkins   Synopsis:
12595305534fSZach Atkins   #include <petscoptions.h>
12605305534fSZach Atkins   PetscErrorCode PetscOptionsIntArray(const char opt[], const char text[], const char man[], PetscInt value[], PetscInt *n, PetscBool *set)
12615305534fSZach Atkins 
12625305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
12635305534fSZach Atkins 
12645305534fSZach Atkins   Input Parameters:
12655305534fSZach Atkins + opt  - the option one is seeking
12665305534fSZach Atkins . text - short string describing option
12675305534fSZach Atkins . man  - manual page for option
12685305534fSZach Atkins - n    - maximum number of values
12695305534fSZach Atkins 
12705305534fSZach Atkins   Output Parameters:
12715305534fSZach Atkins + value - location to copy values
12725305534fSZach Atkins . n     - actual number of values found
12735305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
12745305534fSZach Atkins 
12755305534fSZach Atkins   Level: beginner
12765305534fSZach Atkins 
12775305534fSZach Atkins   Notes:
12785305534fSZach Atkins   The array can be passed as
12795305534fSZach Atkins +   a comma separated list -                                  0,1,2,3,4,5,6,7
12805305534fSZach Atkins .   a range (start\-end+1) -                                  0-8
12815305534fSZach Atkins .   a range with given increment (start\-end+1:inc) -         0-7:2
12825305534fSZach Atkins -   a combination of values and ranges separated by commas -  0,1-8,8-15:2
12835305534fSZach Atkins 
12845305534fSZach Atkins   There must be no intervening spaces between the values.
12855305534fSZach Atkins 
1286af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
12875305534fSZach Atkins 
12885305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
12895305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
12905305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
12915305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
12925305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
12935305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
12945305534fSZach Atkins M*/
1295c857d48bSBarry Smith   #define PetscOptionsIntArray(opt, text, man, value, n, set)                       PetscOptionsIntArray_Private(PetscOptionsObject, opt, text, man, value, n, set)
12965305534fSZach Atkins 
12975305534fSZach Atkins /*MC
12985305534fSZach Atkins   PetscOptionsStringArray - Gets an array of string values for a particular
12995305534fSZach Atkins   option in the database. The values must be separated with commas with
13005305534fSZach Atkins   no intervening spaces.
13015305534fSZach Atkins 
13025305534fSZach Atkins   Synopsis:
13035305534fSZach Atkins   #include <petscoptions.h>
13045305534fSZach Atkins   PetscErrorCode PetscOptionsStringArray(const char opt[], const char text[], const char man[], char *value[], PetscInt *nmax, PetscBool  *set)
13055305534fSZach Atkins 
13065305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`; No Fortran Support
13075305534fSZach Atkins 
13085305534fSZach Atkins   Input Parameters:
13095305534fSZach Atkins + opt  - the option one is seeking
13105305534fSZach Atkins . text - short string describing option
13115305534fSZach Atkins . man  - manual page for option
1312c857d48bSBarry Smith - n    - maximum number of strings
13135305534fSZach Atkins 
13145305534fSZach Atkins   Output Parameters:
13155305534fSZach Atkins + value - location to copy strings
1316c857d48bSBarry Smith . n     - actual number of strings found
13175305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
13185305534fSZach Atkins 
13195305534fSZach Atkins   Level: beginner
13205305534fSZach Atkins 
13215305534fSZach Atkins   Notes:
13225305534fSZach Atkins   The user should pass in an array of pointers to char, to hold all the
13235305534fSZach Atkins   strings returned by this function.
13245305534fSZach Atkins 
13255305534fSZach Atkins   The user is responsible for deallocating the strings that are
13265305534fSZach Atkins   returned.
13275305534fSZach Atkins 
1328af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
13295305534fSZach Atkins 
13305305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
13315305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
13325305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
13335305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
13345305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
13355305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
13365305534fSZach Atkins M*/
1337c857d48bSBarry Smith   #define PetscOptionsStringArray(opt, text, man, value, n, set)                    PetscOptionsStringArray_Private(PetscOptionsObject, opt, text, man, value, n, set)
13385305534fSZach Atkins 
13395305534fSZach Atkins /*MC
13405305534fSZach Atkins   PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
13415305534fSZach Atkins   option in the database. The values must be separated with commas with
13425305534fSZach Atkins   no intervening spaces.
13435305534fSZach Atkins 
13445305534fSZach Atkins   Synopsis:
13455305534fSZach Atkins   #include <petscoptions.h>
13465305534fSZach Atkins   PetscErrorCode PetscOptionsBoolArray(const char opt[], const char text[], const char man[], PetscBool value[], PetscInt *n, PetscBool *set)
13475305534fSZach Atkins 
13485305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
13495305534fSZach Atkins 
13505305534fSZach Atkins   Input Parameters:
13515305534fSZach Atkins + opt  - the option one is seeking
13525305534fSZach Atkins . text - short string describing option
13535305534fSZach Atkins . man  - manual page for option
13545305534fSZach Atkins - n    - maximum number of values allowed in the value array
13555305534fSZach Atkins 
13565305534fSZach Atkins   Output Parameters:
13575305534fSZach Atkins + value - location to copy values
13585305534fSZach Atkins . n     - actual number of values found
13595305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
13605305534fSZach Atkins 
13615305534fSZach Atkins   Level: beginner
13625305534fSZach Atkins 
13635305534fSZach Atkins   Notes:
13645305534fSZach Atkins   The user should pass in an array of `PetscBool`
13655305534fSZach Atkins 
1366af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
13675305534fSZach Atkins 
13685305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
13695305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
13705305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
13715305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
13725305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
13735305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
13745305534fSZach Atkins M*/
1375c857d48bSBarry Smith   #define PetscOptionsBoolArray(opt, text, man, value, n, set)                      PetscOptionsBoolArray_Private(PetscOptionsObject, opt, text, man, value, n, set)
13765305534fSZach Atkins 
13775305534fSZach Atkins /*MC
13785305534fSZach Atkins   PetscOptionsEnumArray - Gets an array of enum values for a particular
13795305534fSZach Atkins   option in the database.
13805305534fSZach Atkins 
13815305534fSZach Atkins   Synopsis:
13825305534fSZach Atkins   #include <petscoptions.h>
13835305534fSZach Atkins   PetscErrorCode PetscOptionsEnumArray(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum value[], PetscInt *n, PetscBool *set)
13845305534fSZach Atkins 
13855305534fSZach Atkins   Logically Collective on the communicator passed in `PetscOptionsBegin()`
13865305534fSZach Atkins 
13875305534fSZach Atkins   Input Parameters:
13885305534fSZach Atkins + opt  - the option one is seeking
13895305534fSZach Atkins . text - short string describing option
13905305534fSZach Atkins . man  - manual page for option
13915305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
13925305534fSZach Atkins - n    - maximum number of values allowed in the value array
13935305534fSZach Atkins 
13945305534fSZach Atkins   Output Parameters:
13955305534fSZach Atkins + value - location to copy values
13965305534fSZach Atkins . n     - actual number of values found
13975305534fSZach Atkins - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
13985305534fSZach Atkins 
13995305534fSZach Atkins   Level: beginner
14005305534fSZach Atkins 
14015305534fSZach Atkins   Notes:
14025305534fSZach Atkins   The array must be passed as a comma separated list.
14035305534fSZach Atkins 
14045305534fSZach Atkins   There must be no intervening spaces between the values.
14055305534fSZach Atkins 
1406af27ebaaSBarry Smith   Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`
14075305534fSZach Atkins 
14085305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
14095305534fSZach Atkins           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`,
14105305534fSZach Atkins           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
14115305534fSZach Atkins           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
14125305534fSZach Atkins           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
14135305534fSZach Atkins           `PetscOptionsFList()`, `PetscOptionsEList()`
14145305534fSZach Atkins M*/
14151ff8fb82SZach Atkins   #define PetscOptionsEnumArray(opt, text, man, list, value, n, set)                PetscOptionsEnumArray_Private(PetscOptionsObject, opt, text, man, list, value, n, set)
14165305534fSZach Atkins 
14175305534fSZach Atkins /*MC
14185305534fSZach Atkins   PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with `newname`
14195305534fSZach Atkins 
14205305534fSZach Atkins   Prints a deprecation warning, unless an option is supplied to suppress.
14215305534fSZach Atkins 
14225305534fSZach Atkins   Logically Collective
14235305534fSZach Atkins 
14245305534fSZach Atkins   Input Parameters:
14255305534fSZach Atkins + oldname - the old, deprecated option
14265305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed
14275305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9"
14285305534fSZach Atkins - info    - additional information string, or `NULL`.
14295305534fSZach Atkins 
14305305534fSZach Atkins   Options Database Key:
14315305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings
14325305534fSZach Atkins 
14335305534fSZach Atkins   Level: developer
14345305534fSZach Atkins 
14355305534fSZach Atkins   Notes:
14365305534fSZach Atkins   If `newname` is provided then the options database will automatically check the database for `oldname`.
14375305534fSZach Atkins 
14385305534fSZach Atkins   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
14395305534fSZach Atkins   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
14405305534fSZach Atkins   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
14415305534fSZach Atkins 
14425305534fSZach Atkins   Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
1443af27ebaaSBarry Smith   Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information
14445305534fSZach Atkins   If newname is provided, the old option is replaced. Otherwise, it remains in the options database.
14455305534fSZach Atkins   If an option is not replaced, the info argument should be used to advise the user on how to proceed.
14465305534fSZach Atkins   There is a limit on the length of the warning printed, so very long strings provided as info may be truncated.
14475305534fSZach Atkins 
14485305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
14495305534fSZach Atkins M*/
14501ff8fb82SZach Atkins   #define PetscOptionsDeprecated(opt, text, man, info)                              PetscOptionsDeprecated_Private(PetscOptionsObject, opt, text, man, info)
14515305534fSZach Atkins 
14525305534fSZach Atkins /*MC
14535305534fSZach Atkins   PetscOptionsDeprecatedMoObject - mark an option as deprecated in the global PetscOptionsObject, optionally replacing it with `newname`
14545305534fSZach Atkins 
14555305534fSZach Atkins   Prints a deprecation warning, unless an option is supplied to suppress.
14565305534fSZach Atkins 
14575305534fSZach Atkins   Logically Collective
14585305534fSZach Atkins 
14595305534fSZach Atkins   Input Parameters:
14605305534fSZach Atkins + oldname - the old, deprecated option
14615305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed
14625305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9"
14635305534fSZach Atkins - info    - additional information string, or `NULL`.
14645305534fSZach Atkins 
14655305534fSZach Atkins   Options Database Key:
14665305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings
14675305534fSZach Atkins 
14685305534fSZach Atkins   Level: developer
14695305534fSZach Atkins 
14705305534fSZach Atkins   Notes:
14715305534fSZach Atkins   If `newname` is provided then the options database will automatically check the database for `oldname`.
14725305534fSZach Atkins 
14735305534fSZach Atkins   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
14745305534fSZach Atkins   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
14755305534fSZach Atkins   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
14765305534fSZach Atkins 
1477af27ebaaSBarry Smith   Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information
14785305534fSZach Atkins   If newname is provided, the old option is replaced. Otherwise, it remains in the options database.
14795305534fSZach Atkins   If an option is not replaced, the info argument should be used to advise the user on how to proceed.
14805305534fSZach Atkins   There is a limit on the length of the warning printed, so very long strings provided as info may be truncated.
14815305534fSZach Atkins 
14825305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
14835305534fSZach Atkins M*/
14841ff8fb82SZach Atkins   #define PetscOptionsDeprecatedNoObject(opt, text, man, info)                      PetscOptionsDeprecated_Private(NULL, opt, text, man, info)
14855f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1486e55864a3SBarry Smith 
1487ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems, const char[], const char[], const char[], const char *const *, PetscEnum, PetscEnum *, PetscBool *);
1488ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInt_Private(PetscOptionItems, const char[], const char[], const char[], PetscInt, PetscInt *, PetscBool *, PetscInt, PetscInt);
1489ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsMPIInt_Private(PetscOptionItems, const char[], const char[], const char[], PetscMPIInt, PetscMPIInt *, PetscBool *, PetscMPIInt, PetscMPIInt);
1490ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsReal_Private(PetscOptionItems, const char[], const char[], const char[], PetscReal, PetscReal *, PetscBool *, PetscReal, PetscReal);
1491ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems, const char[], const char[], const char[], PetscScalar, PetscScalar *, PetscBool *);
1492ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsName_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool *);
1493ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsString_Private(PetscOptionItems, const char[], const char[], const char[], const char[], char *, size_t, PetscBool *);
1494ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBool_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool, PetscBool *, PetscBool *);
1495*4d81f786SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBool3_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool3, PetscBool3 *, PetscBool *);
1496ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool *);
1497ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool *);
1498ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool *);
1499ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFList_Private(PetscOptionItems, const char[], const char[], const char[], PetscFunctionList, const char[], char[], size_t, PetscBool *);
1500ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEList_Private(PetscOptionItems, const char[], const char[], const char[], const char *const *, PetscInt, const char[], PetscInt *, PetscBool *);
1501ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems, const char[], const char[], const char[], PetscReal[], PetscInt *, PetscBool *);
1502ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems, const char[], const char[], const char[], PetscScalar[], PetscInt *, PetscBool *);
1503ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems, const char[], const char[], const char[], PetscInt[], PetscInt *, PetscBool *);
1504ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems, const char[], const char[], const char[], char *[], PetscInt *, PetscBool *);
1505ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems, const char[], const char[], const char[], PetscBool[], PetscInt *, PetscBool *);
1506ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems, const char[], const char[], const char[], const char *const *, PetscEnum[], PetscInt *, PetscBool *);
1507ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems, const char[], const char[], const char[], const char[]);
1508cffb1e40SBarry Smith 
1509ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject, PetscErrorCode (*)(PetscObject, PetscOptionItems, void *), PetscErrorCode (*)(PetscObject, void *), void *);
1510ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject, PetscOptionItems);
1511447ac60bSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1512447ac60bSBarry Smith 
1513f4bc716fSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeftError(void);
1514