1f26ada1bSBarry Smith /* 237f753daSBarry Smith Routines to determine options set in the options database. 3f26ada1bSBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 5ac09b921SBarry Smith 62c8e378dSBarry Smith #include <petscsys.h> 7c619b03eSJed Brown #include <petscviewertypes.h> 83a3b2205SBarry Smith 9ac09b921SBarry Smith /* SUBMANSEC = Sys */ 10ac09b921SBarry Smith 119355ec05SMatthew G. Knepley typedef enum { 129355ec05SMatthew G. Knepley PETSC_OPT_CODE, 139355ec05SMatthew G. Knepley PETSC_OPT_COMMAND_LINE, 149355ec05SMatthew G. Knepley PETSC_OPT_FILE, 159355ec05SMatthew G. Knepley PETSC_OPT_ENVIRONMENT, 169355ec05SMatthew G. Knepley NUM_PETSC_OPT_SOURCE 179355ec05SMatthew G. Knepley } PetscOptionSource; 189355ec05SMatthew G. Knepley 19c5c1f447SLisandro Dalcin #define PETSC_MAX_OPTION_NAME 512 20c5929fdfSBarry Smith typedef struct _n_PetscOptions *PetscOptions; 21c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsCreate(PetscOptions *); 22b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPush(PetscOptions); 23b4205f0bSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPop(void); 24c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsDestroy(PetscOptions *); 252d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsCreateDefault(void); 262d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsDestroyDefault(void); 27c5929fdfSBarry Smith 282d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsHasHelp(PetscOptions, PetscBool *); 29c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHasName(PetscOptions, const char[], const char[], PetscBool *); 30c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetBool(PetscOptions, const char[], const char[], PetscBool *, PetscBool *); 312d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetInt(PetscOptions, const char[], const char[], PetscInt *, PetscBool *); 322d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnum(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscBool *); 332d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEList(PetscOptions, const char[], const char[], const char *const *, PetscInt, PetscInt *, PetscBool *); 34c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetReal(PetscOptions, const char[], const char[], PetscReal *, PetscBool *); 35c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalar(PetscOptions, const char[], const char[], PetscScalar *, PetscBool *); 362d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetString(PetscOptions, const char[], const char[], char[], size_t, PetscBool *); 372d747510SLisandro Dalcin 382d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetBoolArray(PetscOptions, const char[], const char[], PetscBool[], PetscInt *, PetscBool *); 392d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetEnumArray(PetscOptions, const char[], const char[], const char *const *, PetscEnum *, PetscInt *, PetscBool *); 40c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetIntArray(PetscOptions, const char[], const char[], PetscInt[], PetscInt *, PetscBool *); 41c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetRealArray(PetscOptions, const char[], const char[], PetscReal[], PetscInt *, PetscBool *); 42c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetScalarArray(PetscOptions, const char[], const char[], PetscScalar[], PetscInt *, PetscBool *); 43c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsGetStringArray(PetscOptions, const char[], const char[], char *[], PetscInt *, PetscBool *); 443a3b2205SBarry Smith 452d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsValidKey(const char[], PetscBool *); 46c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetAlias(PetscOptions, const char[], const char[]); 47c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSetValue(PetscOptions, const char[], const char[]); 48c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClearValue(PetscOptions, const char[]); 492d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsFindPair(PetscOptions, const char[], const char[], const char *[], PetscBool *); 503a3b2205SBarry Smith 512d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsGetAll(PetscOptions, char *[]); 52c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsAllUsed(PetscOptions, PetscInt *); 532d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsUsed(PetscOptions, const char[], PetscBool *); 54c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeft(PetscOptions); 551ab23d95SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftGet(PetscOptions, PetscInt *, char ***, char ***); 565b191818SFande Kong PETSC_EXTERN PetscErrorCode PetscOptionsLeftRestore(PetscOptions, PetscInt *, char ***, char ***); 57c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsView(PetscOptions, PetscViewer); 584b0e389bSBarry Smith 592d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsReject(PetscOptions, const char[], const char[], const char[]); 60c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsert(PetscOptions, int *, char ***, const char[]); 61c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertFile(MPI_Comm, PetscOptions, const char[], PetscBool); 625c23ca1cSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertFileYAML(MPI_Comm, PetscOptions, const char[], PetscBool); 63c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInsertString(PetscOptions, const char[]); 64080f0011SToby Isaac PETSC_EXTERN PetscErrorCode PetscOptionsInsertStringYAML(PetscOptions, const char[]); 65d06005cbSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsInsertArgs(PetscOptions, int, char **); 66c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsClear(PetscOptions); 67c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPush(PetscOptions, const char[]); 68c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsPrefixPop(PetscOptions); 695d0dffe5SBarry Smith 70014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsGetenv(MPI_Comm, const char[], char[], size_t, PetscBool *); 712d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToBool(const char[], PetscBool *); 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToInt(const char[], PetscInt *); 73014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsStringToReal(const char[], PetscReal *); 742d747510SLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscOptionsStringToScalar(const char[], PetscScalar *); 752e8a6d31SBarry Smith 769355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*)(const char[], const char[], PetscOptionSource, void *), void *, PetscErrorCode (*)(void **)); 779355ec05SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsMonitorDefault(const char[], const char[], PetscOptionSource, void *); 78081c24baSBoyana Norris 7984761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectSetOptions(PetscObject, PetscOptions); 8084761bfeSJed Brown PETSC_EXTERN PetscErrorCode PetscObjectGetOptions(PetscObject, PetscOptions *); 8184761bfeSJed Brown 82014dd563SJed Brown PETSC_EXTERN PetscBool PetscOptionsPublish; 83e55864a3SBarry Smith 84e55864a3SBarry Smith /* 85e55864a3SBarry Smith See manual page for PetscOptionsBegin() 864416b707SBarry Smith 874416b707SBarry Smith PetscOptionsItem and PetscOptionsItems are a single option (such as ksp_type) and a collection of such single 884416b707SBarry Smith options being handled with a PetscOptionsBegin/End() 894416b707SBarry Smith 90e55864a3SBarry Smith */ 919371c9d4SSatish Balay typedef enum { 929371c9d4SSatish Balay OPTION_INT, 939371c9d4SSatish Balay OPTION_BOOL, 949371c9d4SSatish Balay OPTION_REAL, 959371c9d4SSatish Balay OPTION_FLIST, 969371c9d4SSatish Balay OPTION_STRING, 979371c9d4SSatish Balay OPTION_REAL_ARRAY, 989371c9d4SSatish Balay OPTION_SCALAR_ARRAY, 999371c9d4SSatish Balay OPTION_HEAD, 1009371c9d4SSatish Balay OPTION_INT_ARRAY, 1019371c9d4SSatish Balay OPTION_ELIST, 1029371c9d4SSatish Balay OPTION_BOOL_ARRAY, 1039371c9d4SSatish Balay OPTION_STRING_ARRAY 1049371c9d4SSatish Balay } PetscOptionType; 1059355ec05SMatthew G. Knepley 1064416b707SBarry Smith typedef struct _n_PetscOptionItem *PetscOptionItem; 1074416b707SBarry Smith struct _n_PetscOptionItem { 108e55864a3SBarry Smith char *option; 109e55864a3SBarry Smith char *text; 110e55864a3SBarry Smith void *data; /* used to hold the default value and then any value it is changed to by GUI */ 111e55864a3SBarry Smith PetscFunctionList flist; /* used for available values for PetscOptionsList() */ 112e55864a3SBarry Smith const char *const *list; /* used for available values for PetscOptionsEList() */ 113e55864a3SBarry Smith char nlist; /* number of entries in list */ 114e55864a3SBarry Smith char *man; 115e55864a3SBarry Smith size_t arraylength; /* number of entries in data in the case that it is an array (of PetscInt etc) */ 116e55864a3SBarry Smith PetscBool set; /* the user has changed this value in the GUI */ 117e55864a3SBarry Smith PetscOptionType type; 1184416b707SBarry Smith PetscOptionItem next; 119e55864a3SBarry Smith char *pman; 120e55864a3SBarry Smith void *edata; 121e55864a3SBarry Smith }; 122e55864a3SBarry Smith 1234416b707SBarry Smith typedef struct _p_PetscOptionItems { 124e55864a3SBarry Smith PetscInt count; 1254416b707SBarry Smith PetscOptionItem next; 126e55864a3SBarry Smith char *prefix, *pprefix; 127e55864a3SBarry Smith char *title; 128e55864a3SBarry Smith MPI_Comm comm; 129e55864a3SBarry Smith PetscBool printhelp, changedmethod, alreadyprinted; 130e55864a3SBarry Smith PetscObject object; 131c5929fdfSBarry Smith PetscOptions options; 1324416b707SBarry Smith } PetscOptionItems; 13330de9b25SBarry Smith 1345f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 13594bad497SJacob Faibussowitsch extern PetscOptionItems *PetscOptionsObject; /* declare this so that the PetscOptions stubs work */ 1365f80ce2aSJacob Faibussowitsch PetscErrorCode PetscOptionsBegin(MPI_Comm, const char *, const char *, const char *); 1375f80ce2aSJacob Faibussowitsch PetscErrorCode PetscObjectOptionsBegin(PetscObject); 1385f80ce2aSJacob Faibussowitsch PetscErrorCode PetscOptionsEnd(void); 1395f80ce2aSJacob Faibussowitsch #else 14030de9b25SBarry Smith /*MC 14130de9b25SBarry Smith PetscOptionsBegin - Begins a set of queries on the options database that are related and should be 1421957e957SBarry Smith displayed on the same window of a GUI that allows the user to set the options interactively. Often one should 14316a05f60SBarry Smith use `PetscObjectOptionsBegin()` rather than this call. 14430de9b25SBarry Smith 145f2ba6396SBarry Smith Synopsis: 146aaa7dc30SBarry Smith #include <petscoptions.h> 147f2ba6396SBarry Smith PetscErrorCode PetscOptionsBegin(MPI_Comm comm,const char prefix[],const char title[],const char mansec[]) 14830de9b25SBarry Smith 149fb455bf4SPatrick Sanan Collective 15030de9b25SBarry Smith 15130de9b25SBarry Smith Input Parameters: 15230de9b25SBarry Smith + comm - communicator that shares GUI 15376280437SVaclav Hapla . prefix - options prefix for all options displayed on window (optional) 15430de9b25SBarry Smith . title - short descriptive text, for example "Krylov Solver Options" 15587497f52SBarry Smith - mansec - section of manual pages for options, for example `KSP` (optional) 15630de9b25SBarry Smith 15730de9b25SBarry Smith Level: intermediate 15830de9b25SBarry Smith 15995452b02SPatrick Sanan Notes: 160d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 161d0609cedSBarry Smith 16287497f52SBarry Smith The set of queries needs to be ended by a call to `PetscOptionsEnd()`. 163fb455bf4SPatrick Sanan 16487497f52SBarry Smith One can add subheadings with `PetscOptionsHeadBegin()`. 16530de9b25SBarry Smith 16695452b02SPatrick Sanan Developer Notes: 16716a05f60SBarry Smith `PetscOptionsPublish` is set in `PetscOptionsCheckInitial_Private()` with `-saws_options`. When `PetscOptionsPublish` is set the 16816a05f60SBarry Smith loop between `PetscOptionsBegin()` and `PetscOptionsEnd()` is run THREE times with `PetscOptionsPublishCount` of values -1,0,1. 16916a05f60SBarry Smith Otherwise the loop is run ONCE with a `PetscOptionsPublishCount` of 1. 17016a05f60SBarry Smith + \-1 - `PetscOptionsInt()` etc. just call `PetscOptionsGetInt()` etc. 17116a05f60SBarry Smith . 0 - The GUI objects are created in `PetscOptionsInt()` etc. and displayed in `PetscOptionsEnd()` and the options 17216a05f60SBarry Smith database updated with user changes; `PetscOptionsGetInt()` etc. are also called. 17316a05f60SBarry Smith - 1 - `PetscOptionsInt()` etc. again call `PetscOptionsGetInt()` etc. (possibly getting new values), in addition the help message and 174fb455bf4SPatrick Sanan default values are printed if -help was given. 17516a05f60SBarry Smith When `PetscOptionsObject.changedmethod` is set this causes `PetscOptionsPublishCount` to be reset to -2 (so in the next loop iteration it is -1) 17616a05f60SBarry Smith and the whole process is repeated. This is to handle when, for example, the `KSPType` is changed thus changing the list of 17716a05f60SBarry Smith options available so they need to be redisplayed so the user can change the. Changing `PetscOptionsObjects.changedmethod` is never 178fb455bf4SPatrick Sanan currently set. 179aee2cecaSBarry Smith 1804bb2516aSBarry Smith Fortran Note: 1814bb2516aSBarry Smith Returns ierr error code per PETSc Fortran API 1824bb2516aSBarry Smith 183db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 184db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 185db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 186db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 187c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 188db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 189db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()` 19030de9b25SBarry Smith M*/ 1919371c9d4SSatish Balay #define PetscOptionsBegin(comm, prefix, mess, sec) \ 1929371c9d4SSatish Balay do { \ 1934416b707SBarry Smith PetscOptionItems PetscOptionsObjectBase; \ 1944416b707SBarry Smith PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \ 195d0609cedSBarry Smith PetscCall(PetscMemzero(PetscOptionsObject, sizeof(*PetscOptionsObject))); \ 196e55864a3SBarry Smith for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \ 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBegin_Private(PetscOptionsObject, comm, prefix, mess, sec)) 19830de9b25SBarry Smith 1995fefd1ebSJed Brown /*MC 2005fefd1ebSJed Brown PetscObjectOptionsBegin - Begins a set of queries on the options database that are related and should be 2015fefd1ebSJed Brown displayed on the same window of a GUI that allows the user to set the options interactively. 2025fefd1ebSJed Brown 203f2ba6396SBarry Smith Synopsis: 204aaa7dc30SBarry Smith #include <petscoptions.h> 205f2ba6396SBarry Smith PetscErrorCode PetscObjectOptionsBegin(PetscObject obj) 2065fefd1ebSJed Brown 207c3339decSBarry Smith Collective 2085fefd1ebSJed Brown 2092fe279fdSBarry Smith Input Parameter: 2105fefd1ebSJed Brown . obj - object to set options for 2115fefd1ebSJed Brown 2125fefd1ebSJed Brown Level: intermediate 2135fefd1ebSJed Brown 21495452b02SPatrick Sanan Notes: 215d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 216d0609cedSBarry Smith 21787497f52SBarry Smith Needs to be ended by a call the `PetscOptionsEnd()` 218d0609cedSBarry Smith 21987497f52SBarry Smith Can add subheadings with `PetscOptionsHeadBegin()` 2205fefd1ebSJed Brown 221db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 222db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 223db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 224db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 225c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 226db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 227db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 2285fefd1ebSJed Brown M*/ 2299371c9d4SSatish Balay #define PetscObjectOptionsBegin(obj) \ 2309371c9d4SSatish Balay do { \ 2314416b707SBarry Smith PetscOptionItems PetscOptionsObjectBase; \ 2324416b707SBarry Smith PetscOptionItems *PetscOptionsObject = &PetscOptionsObjectBase; \ 233c5929fdfSBarry Smith PetscOptionsObject->options = ((PetscObject)obj)->options; \ 234e55864a3SBarry Smith for (PetscOptionsObject->count = (PetscOptionsPublish ? -1 : 1); PetscOptionsObject->count < 2; PetscOptionsObject->count++) { \ 235dbbe0bcdSBarry Smith PetscCall(PetscObjectOptionsBegin_Private(obj, PetscOptionsObject)) 2363194b578SJed Brown 23730de9b25SBarry Smith /*MC 23830de9b25SBarry Smith PetscOptionsEnd - Ends a set of queries on the options database that are related and should be 23930de9b25SBarry Smith displayed on the same window of a GUI that allows the user to set the options interactively. 24030de9b25SBarry Smith 241f2ba6396SBarry Smith Synopsis: 242aaa7dc30SBarry Smith #include <petscoptions.h> 243f2ba6396SBarry Smith PetscErrorCode PetscOptionsEnd(void) 24430de9b25SBarry Smith 2457cdbe19fSJose E. Roman Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()` 2467cdbe19fSJose E. Roman 24730de9b25SBarry Smith Level: intermediate 24830de9b25SBarry Smith 24995452b02SPatrick Sanan Notes: 25087497f52SBarry Smith Needs to be preceded by a call to `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` 25130de9b25SBarry Smith 252d0609cedSBarry Smith This is a macro that handles its own error checking, it does not return an error code. 253d0609cedSBarry Smith 2544bb2516aSBarry Smith Fortran Note: 2554bb2516aSBarry Smith Returns ierr error code per PETSc Fortran API 2564bb2516aSBarry Smith 257db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 258db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 259db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 260db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsHeadBegin()`, 261c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 262db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 263db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscObjectOptionsBegin()` 26430de9b25SBarry Smith M*/ 2659371c9d4SSatish Balay #define PetscOptionsEnd() \ 2669371c9d4SSatish Balay PetscCall(PetscOptionsEnd_Private(PetscOptionsObject)); \ 2679371c9d4SSatish Balay } \ 2689371c9d4SSatish Balay } \ 2699371c9d4SSatish Balay while (0) 2705f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 27130de9b25SBarry Smith 2724416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *, MPI_Comm, const char[], const char[], const char[]); 273dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject, PetscOptionItems *); 2744416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *); 275d0609cedSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsHeadBegin(PetscOptionItems *, const char[]); 27630de9b25SBarry Smith 2775f80ce2aSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 2789371c9d4SSatish Balay template <typename... T> 2799371c9d4SSatish Balay void PetscOptionsHeadBegin(T...); 280d0609cedSBarry Smith void PetscOptionsHeadEnd(void); 2819371c9d4SSatish Balay template <typename... T> 2829371c9d4SSatish Balay PetscErrorCode PetscOptionsEnum(T...); 2839371c9d4SSatish Balay template <typename... T> 2849371c9d4SSatish Balay PetscErrorCode PetscOptionsInt(T...); 2859371c9d4SSatish Balay template <typename... T> 2869371c9d4SSatish Balay PetscErrorCode PetscOptionsBoundedInt(T...); 2879371c9d4SSatish Balay template <typename... T> 2889371c9d4SSatish Balay PetscErrorCode PetscOptionsRangeInt(T...); 2899371c9d4SSatish Balay template <typename... T> 2909371c9d4SSatish Balay PetscErrorCode PetscOptionsReal(T...); 2919371c9d4SSatish Balay template <typename... T> 2929371c9d4SSatish Balay PetscErrorCode PetscOptionsScalar(T...); 2939371c9d4SSatish Balay template <typename... T> 2949371c9d4SSatish Balay PetscErrorCode PetscOptionsName(T...); 2959371c9d4SSatish Balay template <typename... T> 2969371c9d4SSatish Balay PetscErrorCode PetscOptionsString(T...); 2979371c9d4SSatish Balay template <typename... T> 2989371c9d4SSatish Balay PetscErrorCode PetscOptionsBool(T...); 2999371c9d4SSatish Balay template <typename... T> 3009371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(T...); 3019371c9d4SSatish Balay template <typename... T> 3029371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroup(T...); 3039371c9d4SSatish Balay template <typename... T> 3049371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(T...); 3059371c9d4SSatish Balay template <typename... T> 3069371c9d4SSatish Balay PetscErrorCode PetscOptionsFList(T...); 3079371c9d4SSatish Balay template <typename... T> 3089371c9d4SSatish Balay PetscErrorCode PetscOptionsEList(T...); 3099371c9d4SSatish Balay template <typename... T> 3109371c9d4SSatish Balay PetscErrorCode PetscOptionsRealArray(T...); 3119371c9d4SSatish Balay template <typename... T> 3129371c9d4SSatish Balay PetscErrorCode PetscOptionsScalarArray(T...); 3139371c9d4SSatish Balay template <typename... T> 3149371c9d4SSatish Balay PetscErrorCode PetscOptionsIntArray(T...); 3159371c9d4SSatish Balay template <typename... T> 3169371c9d4SSatish Balay PetscErrorCode PetscOptionsStringArray(T...); 3179371c9d4SSatish Balay template <typename... T> 3189371c9d4SSatish Balay PetscErrorCode PetscOptionsBoolArray(T...); 3199371c9d4SSatish Balay template <typename... T> 3209371c9d4SSatish Balay PetscErrorCode PetscOptionsEnumArray(T...); 3219371c9d4SSatish Balay template <typename... T> 3229371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecated(T...); 3239371c9d4SSatish Balay template <typename... T> 3249371c9d4SSatish Balay PetscErrorCode PetscOptionsDeprecatedNoObject(T...); 3255f80ce2aSJacob Faibussowitsch #else 32630de9b25SBarry Smith /*MC 327d0609cedSBarry Smith PetscOptionsHeadBegin - Puts a heading before listing any more published options. Used, for example, 32816a05f60SBarry Smith in `KSPSetFromOptions_GMRES()`. 329d0609cedSBarry Smith 33016a05f60SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 331d0609cedSBarry Smith 332d0609cedSBarry Smith Input Parameter: 333d0609cedSBarry Smith . head - the heading text 334d0609cedSBarry Smith 33587497f52SBarry Smith Level: developer 336d0609cedSBarry Smith 337d0609cedSBarry Smith Notes: 338d0609cedSBarry Smith Handles errors directly, hence does not return an error code 339d0609cedSBarry Smith 34016a05f60SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()`, and `PetscOptionsObject` created in `PetscOptionsBegin()` should be the first argument 341d0609cedSBarry Smith 34287497f52SBarry Smith Must be followed by a call to `PetscOptionsHeadEnd()` in the same function. 343d0609cedSBarry Smith 344db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 345db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 346db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 347c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 348db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 349db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 350a54bb2a9SBarry Smith M*/ 3519371c9d4SSatish Balay #define PetscOptionsHeadBegin(PetscOptionsObject, head) \ 3529371c9d4SSatish Balay do { \ 35348a46eb9SPierre Jolivet if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " %s\n", head)); \ 354d0609cedSBarry Smith } while (0) 355d0609cedSBarry Smith 356edd03b47SJacob Faibussowitsch #define PetscOptionsHead(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadBegin()", ) PetscOptionsHeadBegin(__VA_ARGS__) 357d0609cedSBarry Smith 358d0609cedSBarry Smith /*MC 35987497f52SBarry Smith PetscOptionsHeadEnd - Ends a section of options begun with `PetscOptionsHeadBegin()` 36016a05f60SBarry Smith See, for example, `KSPSetFromOptions_GMRES()`. 36130de9b25SBarry Smith 362f2ba6396SBarry Smith Synopsis: 363aaa7dc30SBarry Smith #include <petscoptions.h> 364d0609cedSBarry Smith PetscErrorCode PetscOptionsHeadEnd(void) 36530de9b25SBarry Smith 3667cdbe19fSJose E. Roman Collective on the comm used in `PetscOptionsBegin()` or obj used in `PetscObjectOptionsBegin()` 3677cdbe19fSJose E. Roman 36830de9b25SBarry Smith Level: intermediate 36930de9b25SBarry Smith 37095452b02SPatrick Sanan Notes: 37187497f52SBarry Smith Must be between a `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` and a `PetscOptionsEnd()` 37230de9b25SBarry Smith 37387497f52SBarry Smith Must be preceded by a call to `PetscOptionsHeadBegin()` in the same function. 37430de9b25SBarry Smith 37587497f52SBarry Smith This needs to be used only if the code below `PetscOptionsHeadEnd()` can be run ONLY once. 37616a05f60SBarry Smith See, for example, `PCSetFromOptions_Composite()`. This is a `return(0)` in it for early exit 377b52f573bSBarry Smith from the function. 378b52f573bSBarry Smith 37956752e42SBarry Smith This is only for use with the PETSc options GUI 380b52f573bSBarry Smith 381db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 382db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 383db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 384c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 385db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 386db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 38730de9b25SBarry Smith M*/ 3889371c9d4SSatish Balay #define PetscOptionsHeadEnd() \ 3899371c9d4SSatish Balay do { \ 3903ba16761SJacob Faibussowitsch if (PetscOptionsObject->count != 1) PetscFunctionReturn(PETSC_SUCCESS); \ 3919371c9d4SSatish Balay } while (0) 392d0609cedSBarry Smith 393edd03b47SJacob Faibussowitsch #define PetscOptionsTail(...) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscOptionsHeadEnd()", ) PetscOptionsHeadEnd(__VA_ARGS__) 394186905e3SBarry Smith 3955305534fSZach Atkins /*MC 3965305534fSZach Atkins PetscOptionsEnum - Gets the enum value for a particular option in the database. 3975305534fSZach Atkins 3985305534fSZach Atkins Synopsis: 3995305534fSZach Atkins #include <petscoptions.h> 4005305534fSZach Atkins PetscErrorCode PetscOptionsEnum(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum currentvalue, PetscEnum *value, PetscBool *set) 4015305534fSZach Atkins 4025305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 4035305534fSZach Atkins 4045305534fSZach Atkins Input Parameters: 4055305534fSZach Atkins + opt - option name 4065305534fSZach Atkins . text - short string that describes the option 4075305534fSZach Atkins . man - manual page with additional information on option 4085305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 4095305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 4105305534fSZach Atkins .vb 4115305534fSZach Atkins PetscOptionsEnum(..., obj->value,&object->value,...) or 4125305534fSZach Atkins value = defaultvalue 4139314d9b7SBarry Smith PetscOptionsEnum(..., value,&value,&set); 4149314d9b7SBarry Smith if (set) { 4155305534fSZach Atkins .ve 4165305534fSZach Atkins 4175305534fSZach Atkins Output Parameters: 4185305534fSZach Atkins + value - the value to return 4195305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 4205305534fSZach Atkins 4215305534fSZach Atkins Level: beginner 4225305534fSZach Atkins 4235305534fSZach Atkins Notes: 4245305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 4255305534fSZach Atkins 4269314d9b7SBarry Smith `list` is usually something like `PCASMTypes` or some other predefined list of enum names 4275305534fSZach Atkins 4285305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 4299314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 4305305534fSZach Atkins 4315305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 4325305534fSZach Atkins 4335305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 4345305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 4355305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 4365305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 4375305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 4385305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 4395305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 4405305534fSZach Atkins M*/ 4411ff8fb82SZach Atkins #define PetscOptionsEnum(opt, text, man, list, currentvalue, value, set) PetscOptionsEnum_Private(PetscOptionsObject, opt, text, man, list, currentvalue, value, set) 4425305534fSZach Atkins 4435305534fSZach Atkins /*MC 4445305534fSZach Atkins PetscOptionsInt - Gets the integer value for a particular option in the database. 4455305534fSZach Atkins 4465305534fSZach Atkins Synopsis: 4475305534fSZach Atkins #include <petscoptions.h> 4485305534fSZach Atkins PetscErrorCode PetscOptionsInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set) 4495305534fSZach Atkins 4505305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 4515305534fSZach Atkins 4525305534fSZach Atkins Input Parameters: 4535305534fSZach Atkins + opt - option name 4545305534fSZach Atkins . text - short string that describes the option 4555305534fSZach Atkins . man - manual page with additional information on option 4565305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 4575305534fSZach Atkins .vb 4585305534fSZach Atkins PetscOptionsInt(..., obj->value, &obj->value, ...) or 4595305534fSZach Atkins value = defaultvalue 4609314d9b7SBarry Smith PetscOptionsInt(..., value, &value, &set); 4619314d9b7SBarry Smith if (set) { 4625305534fSZach Atkins .ve 4635305534fSZach Atkins 4645305534fSZach Atkins Output Parameters: 4655305534fSZach Atkins + value - the integer value to return 4665305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 4675305534fSZach Atkins 4685305534fSZach Atkins Level: beginner 4695305534fSZach Atkins 4705305534fSZach Atkins Notes: 4715305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 4729314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 4735305534fSZach Atkins 4745305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 4755305534fSZach Atkins 4765305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 4775305534fSZach Atkins 4785305534fSZach Atkins .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 4795305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 4805305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 4815305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 4825305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 4835305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 484*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 4855305534fSZach Atkins M*/ 4861ff8fb82SZach Atkins #define PetscOptionsInt(opt, text, man, currentvalue, value, set) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MIN_INT, PETSC_MAX_INT) 4875305534fSZach Atkins 4885305534fSZach Atkins /*MC 489*52ce0ab5SPierre Jolivet PetscOptionsBoundedInt - Gets an integer value greater than or equal to a given bound for a particular option in the database. 4905305534fSZach Atkins 4915305534fSZach Atkins Synopsis: 4925305534fSZach Atkins #include <petscoptions.h> 4939314d9b7SBarry Smith PetscErrorCode PetscOptionsBoundedInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt bound) 4945305534fSZach Atkins 4955305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 4965305534fSZach Atkins 4975305534fSZach Atkins Input Parameters: 4985305534fSZach Atkins + opt - option name 4995305534fSZach Atkins . text - short string that describes the option 5005305534fSZach Atkins . man - manual page with additional information on option 5015305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 5025305534fSZach Atkins .vb 503*52ce0ab5SPierre Jolivet PetscOptionsBoundedInt(..., obj->value, &obj->value, ...) 5045305534fSZach Atkins .ve 5055305534fSZach Atkins or 5065305534fSZach Atkins .vb 5075305534fSZach Atkins value = defaultvalue 508*52ce0ab5SPierre Jolivet PetscOptionsBoundedInt(..., value, &value, &set, ...); 5099314d9b7SBarry Smith if (set) { 5105305534fSZach Atkins .ve 511*52ce0ab5SPierre Jolivet - bound - the requested value should be greater than or equal to this bound or an error is generated 5125305534fSZach Atkins 5135305534fSZach Atkins Output Parameters: 5145305534fSZach Atkins + value - the integer value to return 5159314d9b7SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 5165305534fSZach Atkins 5175305534fSZach Atkins Level: beginner 5185305534fSZach Atkins 5195305534fSZach Atkins Notes: 5205305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 5219314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 5225305534fSZach Atkins 5235305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 5245305534fSZach Atkins 525af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 5265305534fSZach Atkins 5275305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 5285305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 5295305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 5305305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 5315305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 5325305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 533*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 5345305534fSZach Atkins M*/ 5351ff8fb82SZach Atkins #define PetscOptionsBoundedInt(opt, text, man, currentvalue, value, set, lb) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_MAX_INT) 5365305534fSZach Atkins 5375305534fSZach Atkins /*MC 5385305534fSZach Atkins PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database. 5395305534fSZach Atkins 5405305534fSZach Atkins Synopsis: 5415305534fSZach Atkins #include <petscoptions.h> 5429314d9b7SBarry Smith PetscErrorCode PetscOptionsRangeInt(const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt lb, PetscInt ub) 5435305534fSZach Atkins 5445305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 5455305534fSZach Atkins 5465305534fSZach Atkins Input Parameters: 5475305534fSZach Atkins + opt - option name 5485305534fSZach Atkins . text - short string that describes the option 5495305534fSZach Atkins . man - manual page with additional information on option 5505305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 5515305534fSZach Atkins .vb 552*52ce0ab5SPierre Jolivet PetscOptionsRangeInt(..., obj->value, &obj->value, ...) 553*52ce0ab5SPierre Jolivet .ve 554*52ce0ab5SPierre Jolivet or 555*52ce0ab5SPierre Jolivet .vb 5565305534fSZach Atkins value = defaultvalue 557*52ce0ab5SPierre Jolivet PetscOptionsRangeInt(..., value, &value, &set, ...); 5589314d9b7SBarry Smith if (set) { 5595305534fSZach Atkins .ve 5605305534fSZach Atkins . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 5615305534fSZach Atkins - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 5625305534fSZach Atkins 5635305534fSZach Atkins Output Parameters: 5645305534fSZach Atkins + value - the integer value to return 5659314d9b7SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 5665305534fSZach Atkins 5675305534fSZach Atkins Level: beginner 5685305534fSZach Atkins 5695305534fSZach Atkins Notes: 5705305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 5719314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 5725305534fSZach Atkins 5735305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 5745305534fSZach Atkins 575af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 5765305534fSZach Atkins 5775305534fSZach Atkins .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 5785305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 5795305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 5805305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 5815305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 5825305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 583*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 5845305534fSZach Atkins M*/ 5851ff8fb82SZach Atkins #define PetscOptionsRangeInt(opt, text, man, currentvalue, value, set, lb, ub) PetscOptionsInt_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub) 5865305534fSZach Atkins 5875305534fSZach Atkins /*MC 588*52ce0ab5SPierre Jolivet PetscOptionsReal - Gets a `PetscReal` value for a particular option in the database. 5895305534fSZach Atkins 5905305534fSZach Atkins Synopsis: 5915305534fSZach Atkins #include <petscoptions.h> 5925305534fSZach Atkins PetscErrorCode PetscOptionsReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set) 5935305534fSZach Atkins 5945305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 5955305534fSZach Atkins 5965305534fSZach Atkins Input Parameters: 5975305534fSZach Atkins + opt - option name 5985305534fSZach Atkins . text - short string that describes the option 5995305534fSZach Atkins . man - manual page with additional information on option 6005305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6015305534fSZach Atkins .vb 6025305534fSZach Atkins PetscOptionsReal(..., obj->value,&obj->value,...) or 6035305534fSZach Atkins value = defaultvalue 6049314d9b7SBarry Smith PetscOptionsReal(..., value,&value,&set); 6059314d9b7SBarry Smith if (set) { 6065305534fSZach Atkins .ve 6075305534fSZach Atkins 6085305534fSZach Atkins Output Parameters: 6095305534fSZach Atkins + value - the value to return 6105305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 6115305534fSZach Atkins 6125305534fSZach Atkins Level: beginner 6135305534fSZach Atkins 6145305534fSZach Atkins Notes: 6155305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 6169314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 6175305534fSZach Atkins 6185305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 6195305534fSZach Atkins 620af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 6215305534fSZach Atkins 6225305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 6235305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 6245305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 6255305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 6265305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 6275305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 628*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedReal()`, `PetscOptionsRangeReal()` 6295305534fSZach Atkins M*/ 630*52ce0ab5SPierre Jolivet #define PetscOptionsReal(opt, text, man, currentvalue, value, set) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, PETSC_MIN_REAL, PETSC_MAX_REAL) 631*52ce0ab5SPierre Jolivet 632*52ce0ab5SPierre Jolivet /*MC 633*52ce0ab5SPierre Jolivet PetscOptionsBoundedReal - Gets a `PetscReal` value greater than or equal to a given bound for a particular option in the database. 634*52ce0ab5SPierre Jolivet 635*52ce0ab5SPierre Jolivet Synopsis: 636*52ce0ab5SPierre Jolivet #include <petscoptions.h> 637*52ce0ab5SPierre Jolivet PetscErrorCode PetscOptionsBoundedReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal bound) 638*52ce0ab5SPierre Jolivet 639*52ce0ab5SPierre Jolivet Logically Collective on the communicator passed in `PetscOptionsBegin()` 640*52ce0ab5SPierre Jolivet 641*52ce0ab5SPierre Jolivet Input Parameters: 642*52ce0ab5SPierre Jolivet + opt - option name 643*52ce0ab5SPierre Jolivet . text - short string that describes the option 644*52ce0ab5SPierre Jolivet . man - manual page with additional information on option 645*52ce0ab5SPierre Jolivet . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 646*52ce0ab5SPierre Jolivet .vb 647*52ce0ab5SPierre Jolivet PetscOptionsBoundedReal(..., obj->value, &obj->value, ...) 648*52ce0ab5SPierre Jolivet .ve 649*52ce0ab5SPierre Jolivet or 650*52ce0ab5SPierre Jolivet .vb 651*52ce0ab5SPierre Jolivet value = defaultvalue 652*52ce0ab5SPierre Jolivet PetscOptionsBoundedReal(..., value, &value, &set, ...); 653*52ce0ab5SPierre Jolivet if (set) { 654*52ce0ab5SPierre Jolivet .ve 655*52ce0ab5SPierre Jolivet - bound - the requested value should be greater than or equal to this bound or an error is generated 656*52ce0ab5SPierre Jolivet 657*52ce0ab5SPierre Jolivet Output Parameters: 658*52ce0ab5SPierre Jolivet + value - the real value to return 659*52ce0ab5SPierre Jolivet - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 660*52ce0ab5SPierre Jolivet 661*52ce0ab5SPierre Jolivet Level: beginner 662*52ce0ab5SPierre Jolivet 663*52ce0ab5SPierre Jolivet Notes: 664*52ce0ab5SPierre Jolivet If the user does not supply the option at all `value` is NOT changed. Thus 665*52ce0ab5SPierre Jolivet you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 666*52ce0ab5SPierre Jolivet 667*52ce0ab5SPierre Jolivet The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 668*52ce0ab5SPierre Jolivet 669*52ce0ab5SPierre Jolivet Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 670*52ce0ab5SPierre Jolivet 671*52ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 672*52ce0ab5SPierre Jolivet `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 673*52ce0ab5SPierre Jolivet `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 674*52ce0ab5SPierre Jolivet `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 675*52ce0ab5SPierre Jolivet `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 676*52ce0ab5SPierre Jolivet `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 677*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsBoundedInt()`, `PetscOptionsRangeReal()` 678*52ce0ab5SPierre Jolivet M*/ 679*52ce0ab5SPierre Jolivet #define PetscOptionsBoundedReal(opt, text, man, currentvalue, value, set, lb) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, PETSC_MAX_REAL) 680*52ce0ab5SPierre Jolivet 681*52ce0ab5SPierre Jolivet /*MC 682*52ce0ab5SPierre Jolivet PetscOptionsRangeReal - Gets a `PetscReal` value within a range of values for a particular option in the database. 683*52ce0ab5SPierre Jolivet 684*52ce0ab5SPierre Jolivet Synopsis: 685*52ce0ab5SPierre Jolivet #include <petscoptions.h> 686*52ce0ab5SPierre Jolivet PetscErrorCode PetscOptionsRangeReal(const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set, PetscReal lb, PetscReal ub) 687*52ce0ab5SPierre Jolivet 688*52ce0ab5SPierre Jolivet Logically Collective on the communicator passed in `PetscOptionsBegin()` 689*52ce0ab5SPierre Jolivet 690*52ce0ab5SPierre Jolivet Input Parameters: 691*52ce0ab5SPierre Jolivet + opt - option name 692*52ce0ab5SPierre Jolivet . text - short string that describes the option 693*52ce0ab5SPierre Jolivet . man - manual page with additional information on option 694*52ce0ab5SPierre Jolivet . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 695*52ce0ab5SPierre Jolivet .vb 696*52ce0ab5SPierre Jolivet PetscOptionsRangeReal(..., obj->value, &obj->value, ...) 697*52ce0ab5SPierre Jolivet .ve 698*52ce0ab5SPierre Jolivet or 699*52ce0ab5SPierre Jolivet .vb 700*52ce0ab5SPierre Jolivet value = defaultvalue 701*52ce0ab5SPierre Jolivet PetscOptionsRangeReal(..., value, &value, &set, ...); 702*52ce0ab5SPierre Jolivet if (set) { 703*52ce0ab5SPierre Jolivet .ve 704*52ce0ab5SPierre Jolivet . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 705*52ce0ab5SPierre Jolivet - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 706*52ce0ab5SPierre Jolivet 707*52ce0ab5SPierre Jolivet Output Parameters: 708*52ce0ab5SPierre Jolivet + value - the value to return 709*52ce0ab5SPierre Jolivet - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 710*52ce0ab5SPierre Jolivet 711*52ce0ab5SPierre Jolivet Level: beginner 712*52ce0ab5SPierre Jolivet 713*52ce0ab5SPierre Jolivet Notes: 714*52ce0ab5SPierre Jolivet If the user does not supply the option at all `value` is NOT changed. Thus 715*52ce0ab5SPierre Jolivet you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 716*52ce0ab5SPierre Jolivet 717*52ce0ab5SPierre Jolivet The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 718*52ce0ab5SPierre Jolivet 719*52ce0ab5SPierre Jolivet Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 720*52ce0ab5SPierre Jolivet 721*52ce0ab5SPierre Jolivet .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 722*52ce0ab5SPierre Jolivet `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 723*52ce0ab5SPierre Jolivet `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 724*52ce0ab5SPierre Jolivet `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 725*52ce0ab5SPierre Jolivet `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 726*52ce0ab5SPierre Jolivet `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 727*52ce0ab5SPierre Jolivet `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRangeInt()`, `PetscOptionsBoundedReal()` 728*52ce0ab5SPierre Jolivet M*/ 729*52ce0ab5SPierre Jolivet #define PetscOptionsRangeReal(opt, text, man, currentvalue, value, set, lb, ub) PetscOptionsReal_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set, lb, ub) 7305305534fSZach Atkins 7315305534fSZach Atkins /*MC 7325305534fSZach Atkins PetscOptionsScalar - Gets the `PetscScalar` value for a particular option in the database. 7335305534fSZach Atkins 7345305534fSZach Atkins Synopsis: 7355305534fSZach Atkins #include <petscoptions.h> 7365305534fSZach Atkins PetscErrorCode PetscOptionsScalar(const char opt[], const char text[], const char man[], PetscScalar currentvalue, PetscScalar *value, PetscBool *set) 7375305534fSZach Atkins 7385305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 7395305534fSZach Atkins 7405305534fSZach Atkins Input Parameters: 7415305534fSZach Atkins + opt - option name 7425305534fSZach Atkins . text - short string that describes the option 7435305534fSZach Atkins . man - manual page with additional information on option 7445305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 7455305534fSZach Atkins .vb 7465305534fSZach Atkins PetscOptionsScalar(..., obj->value,&obj->value,...) or 7475305534fSZach Atkins value = defaultvalue 7489314d9b7SBarry Smith PetscOptionsScalar(..., value,&value,&set); 7499314d9b7SBarry Smith if (set) { 7505305534fSZach Atkins .ve 7515305534fSZach Atkins 7525305534fSZach Atkins Output Parameters: 7535305534fSZach Atkins + value - the value to return 7545305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 7555305534fSZach Atkins 7565305534fSZach Atkins Level: beginner 7575305534fSZach Atkins 7585305534fSZach Atkins Notes: 7595305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 7609314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 7615305534fSZach Atkins 7625305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 7635305534fSZach Atkins 764af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 7655305534fSZach Atkins 7665305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 7675305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 7685305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 7695305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 7705305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 7715305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 7725305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 7735305534fSZach Atkins M*/ 7741ff8fb82SZach Atkins #define PetscOptionsScalar(opt, text, man, currentvalue, value, set) PetscOptionsScalar_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 7755305534fSZach Atkins 7765305534fSZach Atkins /*MC 7775305534fSZach 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 7785305534fSZach Atkins its value is set to false. 7795305534fSZach Atkins 7805305534fSZach Atkins Synopsis: 7815305534fSZach Atkins #include <petscoptions.h> 7829314d9b7SBarry Smith PetscErrorCode PetscOptionsName(const char opt[], const char text[], const char man[], PetscBool *set) 7835305534fSZach Atkins 7845305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 7855305534fSZach Atkins 7865305534fSZach Atkins Input Parameters: 7875305534fSZach Atkins + opt - option name 7885305534fSZach Atkins . text - short string that describes the option 7895305534fSZach Atkins - man - manual page with additional information on option 7905305534fSZach Atkins 7915305534fSZach Atkins Output Parameter: 7929314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 7935305534fSZach Atkins 7945305534fSZach Atkins Level: beginner 7955305534fSZach Atkins 7965305534fSZach Atkins Note: 797af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 7985305534fSZach Atkins 7995305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8005305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8015305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 8025305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8035305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8045305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8055305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8065305534fSZach Atkins M*/ 8079314d9b7SBarry Smith #define PetscOptionsName(opt, text, man, set) PetscOptionsName_Private(PetscOptionsObject, opt, text, man, set) 8085305534fSZach Atkins 8095305534fSZach Atkins /*MC 8105305534fSZach Atkins PetscOptionsString - Gets the string value for a particular option in the database. 8115305534fSZach Atkins 8125305534fSZach Atkins Synopsis: 8135305534fSZach Atkins #include <petscoptions.h> 8145305534fSZach Atkins PetscErrorCode PetscOptionsString(const char opt[], const char text[], const char man[], const char currentvalue[], char value[], size_t len, PetscBool *set) 8155305534fSZach Atkins 8165305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 8175305534fSZach Atkins 8185305534fSZach Atkins Input Parameters: 8195305534fSZach Atkins + opt - option name 8205305534fSZach Atkins . text - short string that describes the option 8215305534fSZach Atkins . man - manual page with additional information on option 8225305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 8235305534fSZach Atkins - len - length of the result string including null terminator 8245305534fSZach Atkins 8255305534fSZach Atkins Output Parameters: 8265305534fSZach Atkins + value - the value to return 8275305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 8285305534fSZach Atkins 8295305534fSZach Atkins Level: beginner 8305305534fSZach Atkins 8315305534fSZach Atkins Notes: 832af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8335305534fSZach Atkins 8349314d9b7SBarry Smith If the user provided no string (for example `-optionname` `-someotheroption`) `set` is set to `PETSC_TRUE` (and the string is filled with nulls). 8355305534fSZach Atkins 8365305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 8379314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that `set` is `PETSC_TRUE`. 8385305534fSZach Atkins 8395305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 8405305534fSZach Atkins 8415305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8425305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8435305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 8445305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8455305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8465305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8475305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8485305534fSZach Atkins M*/ 8491ff8fb82SZach Atkins #define PetscOptionsString(opt, text, man, currentvalue, value, len, set) PetscOptionsString_Private(PetscOptionsObject, opt, text, man, currentvalue, value, len, set) 8505305534fSZach Atkins 8515305534fSZach Atkins /*MC 8525305534fSZach Atkins PetscOptionsBool - Determines if a particular option is in the database with a true or false 8535305534fSZach Atkins 8545305534fSZach Atkins Synopsis: 8555305534fSZach Atkins #include <petscoptions.h> 8565305534fSZach Atkins PetscErrorCode PetscOptionsBool(const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool *flg, PetscBool *set) 8575305534fSZach Atkins 8585305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 8595305534fSZach Atkins 8605305534fSZach Atkins Input Parameters: 8615305534fSZach Atkins + opt - option name 8625305534fSZach Atkins . text - short string that describes the option 8635305534fSZach Atkins . man - manual page with additional information on option 8645305534fSZach Atkins - currentvalue - the current value 8655305534fSZach Atkins 8665305534fSZach Atkins Output Parameters: 8675305534fSZach Atkins + flg - `PETSC_TRUE` or `PETSC_FALSE` 8685305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 8695305534fSZach Atkins 8705305534fSZach Atkins Level: beginner 8715305534fSZach Atkins 8725305534fSZach Atkins Notes: 8735305534fSZach Atkins TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE` 8745305534fSZach Atkins FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE` 8755305534fSZach Atkins 8769314d9b7SBarry 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` 8775305534fSZach Atkins is equivalent to `-requested_bool true` 8785305534fSZach Atkins 8795305534fSZach Atkins If the user does not supply the option at all `flg` is NOT changed. Thus 8809314d9b7SBarry Smith you should ALWAYS initialize the `flg` variable if you access it without first checking that the `set` flag is `PETSC_TRUE`. 8815305534fSZach Atkins 8825305534fSZach Atkins Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8835305534fSZach Atkins 8845305534fSZach Atkins .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 8855305534fSZach Atkins `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 8865305534fSZach Atkins `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 8875305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 8885305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 8895305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 8905305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 8915305534fSZach Atkins M*/ 8921ff8fb82SZach Atkins #define PetscOptionsBool(opt, text, man, currentvalue, value, set) PetscOptionsBool_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 8935305534fSZach Atkins 8945305534fSZach Atkins /*MC 8955305534fSZach Atkins PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 8965305534fSZach Atkins which at most a single value can be true. 8975305534fSZach Atkins 8985305534fSZach Atkins Synopsis: 8995305534fSZach Atkins #include <petscoptions.h> 9009314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[], const char text[], const char man[], PetscBool *set) 9015305534fSZach Atkins 9025305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9035305534fSZach Atkins 9045305534fSZach Atkins Input Parameters: 9055305534fSZach Atkins + opt - option name 9065305534fSZach Atkins . text - short string that describes the option 9075305534fSZach Atkins - man - manual page with additional information on option 9085305534fSZach Atkins 9095305534fSZach Atkins Output Parameter: 9109314d9b7SBarry Smith . set - whether that option was set or not 9115305534fSZach Atkins 9125305534fSZach Atkins Level: intermediate 9135305534fSZach Atkins 9145305534fSZach Atkins Notes: 915af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9165305534fSZach Atkins 9175305534fSZach Atkins Must be followed by 0 or more `PetscOptionsBoolGroup()`s and `PetscOptionsBoolGroupEnd()` 9185305534fSZach Atkins 9195305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 9205305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 9215305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 9225305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 9235305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 9245305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 9255305534fSZach Atkins M*/ 9269314d9b7SBarry Smith #define PetscOptionsBoolGroupBegin(opt, text, man, set) PetscOptionsBoolGroupBegin_Private(PetscOptionsObject, opt, text, man, set) 9275305534fSZach Atkins 9285305534fSZach Atkins /*MC 9295305534fSZach Atkins PetscOptionsBoolGroup - One in a series of logical queries on the options database for 9305305534fSZach Atkins which at most a single value can be true. 9315305534fSZach Atkins 9325305534fSZach Atkins Synopsis: 9335305534fSZach Atkins #include <petscoptions.h> 9349314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroup(const char opt[], const char text[], const char man[], PetscBool *set) 9355305534fSZach Atkins 9365305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9375305534fSZach Atkins 9385305534fSZach Atkins Input Parameters: 9395305534fSZach Atkins + opt - option name 9405305534fSZach Atkins . text - short string that describes the option 9415305534fSZach Atkins - man - manual page with additional information on option 9425305534fSZach Atkins 9435305534fSZach Atkins Output Parameter: 9449314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 9455305534fSZach Atkins 9465305534fSZach Atkins Level: intermediate 9475305534fSZach Atkins 9485305534fSZach Atkins Notes: 949af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9505305534fSZach Atkins 9515305534fSZach Atkins Must follow a `PetscOptionsBoolGroupBegin()` and preceded a `PetscOptionsBoolGroupEnd()` 9525305534fSZach Atkins 9535305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 9545305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 9555305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 9565305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 9575305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 9585305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 9595305534fSZach Atkins M*/ 9609314d9b7SBarry Smith #define PetscOptionsBoolGroup(opt, text, man, set) PetscOptionsBoolGroup_Private(PetscOptionsObject, opt, text, man, set) 9615305534fSZach Atkins 9625305534fSZach Atkins /*MC 9635305534fSZach Atkins PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 9645305534fSZach Atkins which at most a single value can be true. 9655305534fSZach Atkins 9665305534fSZach Atkins Synopsis: 9675305534fSZach Atkins #include <petscoptions.h> 9689314d9b7SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[], const char text[], const char man[], PetscBool *set) 9695305534fSZach Atkins 9705305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 9715305534fSZach Atkins 9725305534fSZach Atkins Input Parameters: 9735305534fSZach Atkins + opt - option name 9745305534fSZach Atkins . text - short string that describes the option 9755305534fSZach Atkins - man - manual page with additional information on option 9765305534fSZach Atkins 9775305534fSZach Atkins Output Parameter: 9789314d9b7SBarry Smith . set - `PETSC_TRUE` if found, else `PETSC_FALSE` 9795305534fSZach Atkins 9805305534fSZach Atkins Level: intermediate 9815305534fSZach Atkins 9825305534fSZach Atkins Notes: 983af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 9845305534fSZach Atkins 9855305534fSZach Atkins Must follow a `PetscOptionsBoolGroupBegin()` 9865305534fSZach Atkins 9875305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 9885305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 9895305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 9905305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 9915305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 9925305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 9935305534fSZach Atkins M*/ 9949314d9b7SBarry Smith #define PetscOptionsBoolGroupEnd(opt, text, man, set) PetscOptionsBoolGroupEnd_Private(PetscOptionsObject, opt, text, man, set) 9955305534fSZach Atkins 9965305534fSZach Atkins /*MC 9975305534fSZach Atkins PetscOptionsFList - Puts a list of option values that a single one may be selected from 9985305534fSZach Atkins 9995305534fSZach Atkins Synopsis: 10005305534fSZach Atkins #include <petscoptions.h> 10015305534fSZach Atkins PetscErrorCode PetscOptionsFList(const char opt[], const char ltext[], const char man[], PetscFunctionList list, const char currentvalue[], char value[], size_t len, PetscBool *set) 10025305534fSZach Atkins 10035305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 10045305534fSZach Atkins 10055305534fSZach Atkins Input Parameters: 10065305534fSZach Atkins + opt - option name 10075305534fSZach Atkins . ltext - short string that describes the option 10085305534fSZach Atkins . man - manual page with additional information on option 10095305534fSZach Atkins . list - the possible choices 10105305534fSZach Atkins . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10115305534fSZach Atkins .vb 10129314d9b7SBarry Smith PetscOptionsFlist(..., obj->value,value,len,&set); 10139314d9b7SBarry Smith if (set) { 10145305534fSZach Atkins .ve 10155305534fSZach Atkins - len - the length of the character array value 10165305534fSZach Atkins 10175305534fSZach Atkins Output Parameters: 10185305534fSZach Atkins + value - the value to return 10195305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 10205305534fSZach Atkins 10215305534fSZach Atkins Level: intermediate 10225305534fSZach Atkins 10235305534fSZach Atkins Notes: 1024af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 10255305534fSZach Atkins 10265305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 10279314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`. 10285305534fSZach Atkins 10295305534fSZach Atkins The `currentvalue` passed into this routine does not get transferred to the output `value` variable automatically. 10305305534fSZach Atkins 10315305534fSZach Atkins See `PetscOptionsEList()` for when the choices are given in a string array 10325305534fSZach Atkins 10335305534fSZach Atkins To get a listing of all currently specified options, 10345305534fSZach Atkins see `PetscOptionsView()` or `PetscOptionsGetAll()` 10355305534fSZach Atkins 1036af27ebaaSBarry Smith Developer Note: 10375305534fSZach Atkins This cannot check for invalid selection because of things like `MATAIJ` that are not included in the list 10385305534fSZach Atkins 10395305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 10405305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 10415305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 10425305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 10435305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 10445305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 10455305534fSZach Atkins M*/ 10461ff8fb82SZach Atkins #define PetscOptionsFList(opt, ltext, man, list, currentvalue, value, len, set) PetscOptionsFList_Private(PetscOptionsObject, opt, ltext, man, list, currentvalue, value, len, set) 10475305534fSZach Atkins 10485305534fSZach Atkins /*MC 10495305534fSZach Atkins PetscOptionsEList - Puts a list of option values that a single one may be selected from 10505305534fSZach Atkins 10515305534fSZach Atkins Synopsis: 10525305534fSZach Atkins #include <petscoptions.h> 10535305534fSZach 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) 10545305534fSZach Atkins 10555305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 10565305534fSZach Atkins 10575305534fSZach Atkins Input Parameters: 10585305534fSZach Atkins + opt - option name 10595305534fSZach Atkins . ltext - short string that describes the option 10605305534fSZach Atkins . man - manual page with additional information on option 10615305534fSZach Atkins . list - the possible choices (one of these must be selected, anything else is invalid) 10625305534fSZach Atkins . ntext - number of choices 10635305534fSZach Atkins - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10645305534fSZach Atkins .vb 10659314d9b7SBarry Smith PetscOptionsEList(..., obj->value,&value,&set); 10669314d9b7SBarry Smith .ve if (set) { 10675305534fSZach Atkins 10685305534fSZach Atkins Output Parameters: 10695305534fSZach Atkins + value - the index of the value to return 10705305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 10715305534fSZach Atkins 10725305534fSZach Atkins Level: intermediate 10735305534fSZach Atkins 10745305534fSZach Atkins Notes: 1075af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 10765305534fSZach Atkins 10775305534fSZach Atkins If the user does not supply the option at all `value` is NOT changed. Thus 10789314d9b7SBarry Smith you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is `PETSC_TRUE`. 10795305534fSZach Atkins 10805305534fSZach Atkins See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList()` 10815305534fSZach Atkins 10825305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 10835305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 10845305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 10855305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 10865305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 10875305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEnum()` 10885305534fSZach Atkins M*/ 10891ff8fb82SZach Atkins #define PetscOptionsEList(opt, ltext, man, list, ntext, currentvalue, value, set) PetscOptionsEList_Private(PetscOptionsObject, opt, ltext, man, list, ntext, currentvalue, value, set) 10905305534fSZach Atkins 10915305534fSZach Atkins /*MC 10925305534fSZach Atkins PetscOptionsRealArray - Gets an array of double values for a particular 10935305534fSZach Atkins option in the database. The values must be separated with commas with 10945305534fSZach Atkins no intervening spaces. 10955305534fSZach Atkins 10965305534fSZach Atkins Synopsis: 10975305534fSZach Atkins #include <petscoptions.h> 10985305534fSZach Atkins PetscErrorCode PetscOptionsRealArray(const char opt[], const char text[], const char man[], PetscReal value[], PetscInt *n, PetscBool *set) 10995305534fSZach Atkins 11005305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11015305534fSZach Atkins 11025305534fSZach Atkins Input Parameters: 11035305534fSZach Atkins + opt - the option one is seeking 11045305534fSZach Atkins . text - short string describing option 11055305534fSZach Atkins . man - manual page for option 11065305534fSZach Atkins - n - maximum number of values that value has room for 11075305534fSZach Atkins 11085305534fSZach Atkins Output Parameters: 11095305534fSZach Atkins + value - location to copy values 11105305534fSZach Atkins . n - actual number of values found 11115305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11125305534fSZach Atkins 11135305534fSZach Atkins Level: beginner 11145305534fSZach Atkins 11155305534fSZach Atkins Note: 1116af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 11175305534fSZach Atkins 11185305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 11195305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 11205305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 11215305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 11225305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 11235305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 11245305534fSZach Atkins M*/ 11251ff8fb82SZach Atkins #define PetscOptionsRealArray(opt, text, man, currentvalue, value, set) PetscOptionsRealArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 11265305534fSZach Atkins 11275305534fSZach Atkins /*MC 11285305534fSZach Atkins PetscOptionsScalarArray - Gets an array of `PetscScalar` values for a particular 11295305534fSZach Atkins option in the database. The values must be separated with commas with 11305305534fSZach Atkins no intervening spaces. 11315305534fSZach Atkins 11325305534fSZach Atkins Synopsis: 11335305534fSZach Atkins #include <petscoptions.h> 11345305534fSZach Atkins PetscErrorCode PetscOptionsScalarArray(const char opt[], const char text[], const char man[], PetscScalar value[], PetscInt *n, PetscBool *set) 11355305534fSZach Atkins 11365305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11375305534fSZach Atkins 11385305534fSZach Atkins Input Parameters: 11395305534fSZach Atkins + opt - the option one is seeking 11405305534fSZach Atkins . text - short string describing option 11415305534fSZach Atkins . man - manual page for option 11425305534fSZach Atkins - n - maximum number of values allowed in the value array 11435305534fSZach Atkins 11445305534fSZach Atkins Output Parameters: 11455305534fSZach Atkins + value - location to copy values 11465305534fSZach Atkins . n - actual number of values found 11475305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11485305534fSZach Atkins 11495305534fSZach Atkins Level: beginner 11505305534fSZach Atkins 11515305534fSZach Atkins Note: 1152af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 11535305534fSZach Atkins 11545305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 11555305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 11565305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 11575305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 11585305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 11595305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 11605305534fSZach Atkins M*/ 11611ff8fb82SZach Atkins #define PetscOptionsScalarArray(opt, text, man, currentvalue, value, set) PetscOptionsScalarArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 11625305534fSZach Atkins 11635305534fSZach Atkins /*MC 11645305534fSZach Atkins PetscOptionsIntArray - Gets an array of integers for a particular 11655305534fSZach Atkins option in the database. 11665305534fSZach Atkins 11675305534fSZach Atkins Synopsis: 11685305534fSZach Atkins #include <petscoptions.h> 11695305534fSZach Atkins PetscErrorCode PetscOptionsIntArray(const char opt[], const char text[], const char man[], PetscInt value[], PetscInt *n, PetscBool *set) 11705305534fSZach Atkins 11715305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 11725305534fSZach Atkins 11735305534fSZach Atkins Input Parameters: 11745305534fSZach Atkins + opt - the option one is seeking 11755305534fSZach Atkins . text - short string describing option 11765305534fSZach Atkins . man - manual page for option 11775305534fSZach Atkins - n - maximum number of values 11785305534fSZach Atkins 11795305534fSZach Atkins Output Parameters: 11805305534fSZach Atkins + value - location to copy values 11815305534fSZach Atkins . n - actual number of values found 11825305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 11835305534fSZach Atkins 11845305534fSZach Atkins Level: beginner 11855305534fSZach Atkins 11865305534fSZach Atkins Notes: 11875305534fSZach Atkins The array can be passed as 11885305534fSZach Atkins + a comma separated list - 0,1,2,3,4,5,6,7 11895305534fSZach Atkins . a range (start\-end+1) - 0-8 11905305534fSZach Atkins . a range with given increment (start\-end+1:inc) - 0-7:2 11915305534fSZach Atkins - a combination of values and ranges separated by commas - 0,1-8,8-15:2 11925305534fSZach Atkins 11935305534fSZach Atkins There must be no intervening spaces between the values. 11945305534fSZach Atkins 1195af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 11965305534fSZach Atkins 11975305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 11985305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 11995305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12005305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12015305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12025305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12035305534fSZach Atkins M*/ 12041ff8fb82SZach Atkins #define PetscOptionsIntArray(opt, text, man, currentvalue, value, set) PetscOptionsIntArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12055305534fSZach Atkins 12065305534fSZach Atkins /*MC 12075305534fSZach Atkins PetscOptionsStringArray - Gets an array of string values for a particular 12085305534fSZach Atkins option in the database. The values must be separated with commas with 12095305534fSZach Atkins no intervening spaces. 12105305534fSZach Atkins 12115305534fSZach Atkins Synopsis: 12125305534fSZach Atkins #include <petscoptions.h> 12135305534fSZach Atkins PetscErrorCode PetscOptionsStringArray(const char opt[], const char text[], const char man[], char *value[], PetscInt *nmax, PetscBool *set) 12145305534fSZach Atkins 12155305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()`; No Fortran Support 12165305534fSZach Atkins 12175305534fSZach Atkins Input Parameters: 12185305534fSZach Atkins + opt - the option one is seeking 12195305534fSZach Atkins . text - short string describing option 12205305534fSZach Atkins . man - manual page for option 12215305534fSZach Atkins - nmax - maximum number of strings 12225305534fSZach Atkins 12235305534fSZach Atkins Output Parameters: 12245305534fSZach Atkins + value - location to copy strings 12255305534fSZach Atkins . nmax - actual number of strings found 12265305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 12275305534fSZach Atkins 12285305534fSZach Atkins Level: beginner 12295305534fSZach Atkins 12305305534fSZach Atkins Notes: 12315305534fSZach Atkins The user should pass in an array of pointers to char, to hold all the 12325305534fSZach Atkins strings returned by this function. 12335305534fSZach Atkins 12345305534fSZach Atkins The user is responsible for deallocating the strings that are 12355305534fSZach Atkins returned. 12365305534fSZach Atkins 1237af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 12385305534fSZach Atkins 12395305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 12405305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 12415305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12425305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12435305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12445305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12455305534fSZach Atkins M*/ 12461ff8fb82SZach Atkins #define PetscOptionsStringArray(opt, text, man, currentvalue, value, set) PetscOptionsStringArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12475305534fSZach Atkins 12485305534fSZach Atkins /*MC 12495305534fSZach Atkins PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 12505305534fSZach Atkins option in the database. The values must be separated with commas with 12515305534fSZach Atkins no intervening spaces. 12525305534fSZach Atkins 12535305534fSZach Atkins Synopsis: 12545305534fSZach Atkins #include <petscoptions.h> 12555305534fSZach Atkins PetscErrorCode PetscOptionsBoolArray(const char opt[], const char text[], const char man[], PetscBool value[], PetscInt *n, PetscBool *set) 12565305534fSZach Atkins 12575305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 12585305534fSZach Atkins 12595305534fSZach Atkins Input Parameters: 12605305534fSZach Atkins + opt - the option one is seeking 12615305534fSZach Atkins . text - short string describing option 12625305534fSZach Atkins . man - manual page for option 12635305534fSZach Atkins - n - maximum number of values allowed in the value array 12645305534fSZach Atkins 12655305534fSZach Atkins Output Parameters: 12665305534fSZach Atkins + value - location to copy values 12675305534fSZach Atkins . n - actual number of values found 12685305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 12695305534fSZach Atkins 12705305534fSZach Atkins Level: beginner 12715305534fSZach Atkins 12725305534fSZach Atkins Notes: 12735305534fSZach Atkins The user should pass in an array of `PetscBool` 12745305534fSZach Atkins 1275af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 12765305534fSZach Atkins 12775305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 12785305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 12795305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 12805305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 12815305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 12825305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 12835305534fSZach Atkins M*/ 12841ff8fb82SZach Atkins #define PetscOptionsBoolArray(opt, text, man, currentvalue, value, set) PetscOptionsBoolArray_Private(PetscOptionsObject, opt, text, man, currentvalue, value, set) 12855305534fSZach Atkins 12865305534fSZach Atkins /*MC 12875305534fSZach Atkins PetscOptionsEnumArray - Gets an array of enum values for a particular 12885305534fSZach Atkins option in the database. 12895305534fSZach Atkins 12905305534fSZach Atkins Synopsis: 12915305534fSZach Atkins #include <petscoptions.h> 12925305534fSZach Atkins PetscErrorCode PetscOptionsEnumArray(const char opt[], const char text[], const char man[], const char *const *list, PetscEnum value[], PetscInt *n, PetscBool *set) 12935305534fSZach Atkins 12945305534fSZach Atkins Logically Collective on the communicator passed in `PetscOptionsBegin()` 12955305534fSZach Atkins 12965305534fSZach Atkins Input Parameters: 12975305534fSZach Atkins + opt - the option one is seeking 12985305534fSZach Atkins . text - short string describing option 12995305534fSZach Atkins . man - manual page for option 13005305534fSZach Atkins . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 13015305534fSZach Atkins - n - maximum number of values allowed in the value array 13025305534fSZach Atkins 13035305534fSZach Atkins Output Parameters: 13045305534fSZach Atkins + value - location to copy values 13055305534fSZach Atkins . n - actual number of values found 13065305534fSZach Atkins - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 13075305534fSZach Atkins 13085305534fSZach Atkins Level: beginner 13095305534fSZach Atkins 13105305534fSZach Atkins Notes: 13115305534fSZach Atkins The array must be passed as a comma separated list. 13125305534fSZach Atkins 13135305534fSZach Atkins There must be no intervening spaces between the values. 13145305534fSZach Atkins 1315af27ebaaSBarry Smith Must be used between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 13165305534fSZach Atkins 13175305534fSZach Atkins .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 13185305534fSZach Atkins `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 13195305534fSZach Atkins `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 13205305534fSZach Atkins `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 13215305534fSZach Atkins `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 13225305534fSZach Atkins `PetscOptionsFList()`, `PetscOptionsEList()` 13235305534fSZach Atkins M*/ 13241ff8fb82SZach Atkins #define PetscOptionsEnumArray(opt, text, man, list, value, n, set) PetscOptionsEnumArray_Private(PetscOptionsObject, opt, text, man, list, value, n, set) 13255305534fSZach Atkins 13265305534fSZach Atkins /*MC 13275305534fSZach Atkins PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with `newname` 13285305534fSZach Atkins 13295305534fSZach Atkins Prints a deprecation warning, unless an option is supplied to suppress. 13305305534fSZach Atkins 13315305534fSZach Atkins Logically Collective 13325305534fSZach Atkins 13335305534fSZach Atkins Input Parameters: 13345305534fSZach Atkins + oldname - the old, deprecated option 13355305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed 13365305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9" 13375305534fSZach Atkins - info - additional information string, or `NULL`. 13385305534fSZach Atkins 13395305534fSZach Atkins Options Database Key: 13405305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings 13415305534fSZach Atkins 13425305534fSZach Atkins Level: developer 13435305534fSZach Atkins 13445305534fSZach Atkins Notes: 13455305534fSZach Atkins If `newname` is provided then the options database will automatically check the database for `oldname`. 13465305534fSZach Atkins 13475305534fSZach Atkins The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the 13485305534fSZach Atkins new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call. 13495305534fSZach Atkins See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails. 13505305534fSZach Atkins 13515305534fSZach Atkins Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`. 1352af27ebaaSBarry Smith Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information 13535305534fSZach Atkins If newname is provided, the old option is replaced. Otherwise, it remains in the options database. 13545305534fSZach Atkins If an option is not replaced, the info argument should be used to advise the user on how to proceed. 13555305534fSZach Atkins There is a limit on the length of the warning printed, so very long strings provided as info may be truncated. 13565305534fSZach Atkins 13575305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()` 13585305534fSZach Atkins M*/ 13591ff8fb82SZach Atkins #define PetscOptionsDeprecated(opt, text, man, info) PetscOptionsDeprecated_Private(PetscOptionsObject, opt, text, man, info) 13605305534fSZach Atkins 13615305534fSZach Atkins /*MC 13625305534fSZach Atkins PetscOptionsDeprecatedMoObject - mark an option as deprecated in the global PetscOptionsObject, optionally replacing it with `newname` 13635305534fSZach Atkins 13645305534fSZach Atkins Prints a deprecation warning, unless an option is supplied to suppress. 13655305534fSZach Atkins 13665305534fSZach Atkins Logically Collective 13675305534fSZach Atkins 13685305534fSZach Atkins Input Parameters: 13695305534fSZach Atkins + oldname - the old, deprecated option 13705305534fSZach Atkins . newname - the new option, or `NULL` if option is purely removed 13715305534fSZach Atkins . version - a string describing the version of first deprecation, e.g. "3.9" 13725305534fSZach Atkins - info - additional information string, or `NULL`. 13735305534fSZach Atkins 13745305534fSZach Atkins Options Database Key: 13755305534fSZach Atkins . -options_suppress_deprecated_warnings - do not print deprecation warnings 13765305534fSZach Atkins 13775305534fSZach Atkins Level: developer 13785305534fSZach Atkins 13795305534fSZach Atkins Notes: 13805305534fSZach Atkins If `newname` is provided then the options database will automatically check the database for `oldname`. 13815305534fSZach Atkins 13825305534fSZach Atkins The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the 13835305534fSZach Atkins new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call. 13845305534fSZach Atkins See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails. 13855305534fSZach Atkins 1386af27ebaaSBarry Smith Only the process of MPI rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or `PetscObjectOptionsBegin()` prints the information 13875305534fSZach Atkins If newname is provided, the old option is replaced. Otherwise, it remains in the options database. 13885305534fSZach Atkins If an option is not replaced, the info argument should be used to advise the user on how to proceed. 13895305534fSZach Atkins There is a limit on the length of the warning printed, so very long strings provided as info may be truncated. 13905305534fSZach Atkins 13915305534fSZach Atkins .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()` 13925305534fSZach Atkins M*/ 13931ff8fb82SZach Atkins #define PetscOptionsDeprecatedNoObject(opt, text, man, info) PetscOptionsDeprecated_Private(NULL, opt, text, man, info) 13945f80ce2aSJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1395e55864a3SBarry Smith 13964416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum, PetscEnum *, PetscBool *); 13975a856986SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt, PetscInt *, PetscBool *, PetscInt, PetscInt); 1398*52ce0ab5SPierre Jolivet PETSC_EXTERN PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal, PetscReal *, PetscBool *, PetscReal, PetscReal); 13994416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar, PetscScalar *, PetscBool *); 14004416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsName_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14014416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsString_Private(PetscOptionItems *, const char[], const char[], const char[], const char[], char *, size_t, PetscBool *); 14024416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool, PetscBool *, PetscBool *); 14034416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14044416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14054416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool *); 14064416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFList_Private(PetscOptionItems *, const char[], const char[], const char[], PetscFunctionList, const char[], char[], size_t, PetscBool *); 14074416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEList_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscInt, const char[], PetscInt *, PetscBool *); 14084416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscReal[], PetscInt *, PetscBool *); 14094416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscScalar[], PetscInt *, PetscBool *); 14104416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscInt[], PetscInt *, PetscBool *); 14114416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *, const char[], const char[], const char[], char *[], PetscInt *, PetscBool *); 14124416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *, const char[], const char[], const char[], PetscBool[], PetscInt *, PetscBool *); 14134416b707SBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *, const char[], const char[], const char[], const char *const *, PetscEnum[], PetscInt *, PetscBool *); 14149f3a6782SPatrick Sanan PETSC_EXTERN PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *, const char[], const char[], const char[], const char[]); 1415cffb1e40SBarry Smith 1416e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsSAWsDestroy(void); 1417f8d0b74dSMatthew Knepley 1418dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject, PetscErrorCode (*)(PetscObject, PetscOptionItems *, void *), PetscErrorCode (*)(PetscObject, void *), void *); 1419dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject, PetscOptionItems *); 1420447ac60bSBarry Smith PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject); 1421447ac60bSBarry Smith 1422f4bc716fSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsLeftError(void); 1423