xref: /petsc/src/sys/objects/options.c (revision 10c654e6b2f9e2edc1293b3f0e047bdcfaeccab6)
173fca5a0SBarry Smith /* Define Feature test macros to make sure atoll is available (SVr4, POSIX.1-2001, 4.3BSD, C99), not in (C89 and POSIX.1-1996) */
20039db0dSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for atoll() */
3e5ea902fSJed Brown 
4e5c89e4eSSatish Balay /*
53fc1eb6aSBarry Smith    These routines simplify the use of command line, file options, etc., and are used to manipulate the options database.
63fc1eb6aSBarry Smith    This provides the low-level interface, the high level interface is in aoptions.c
7e5c89e4eSSatish Balay 
83fc1eb6aSBarry Smith    Some routines use regular malloc and free because it cannot know  what malloc is requested with the
93fc1eb6aSBarry Smith    options database until it has already processed the input.
10e5c89e4eSSatish Balay */
11e5c89e4eSSatish Balay 
12af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I  "petscsys.h"   I*/
13665c2dedSJed Brown #include <petscviewer.h>
14ad1ac5ecSJed Brown #include <ctype.h>
15e5c89e4eSSatish Balay #if defined(PETSC_HAVE_MALLOC_H)
16e5c89e4eSSatish Balay   #include <malloc.h>
17e5c89e4eSSatish Balay #endif
18ef279fd6SBarry Smith #if defined(PETSC_HAVE_STRINGS_H)
19ef279fd6SBarry Smith   #include <strings.h> /* strcasecmp */
20ef279fd6SBarry Smith #endif
21e5c89e4eSSatish Balay 
222d747510SLisandro Dalcin #if defined(PETSC_HAVE_STRCASECMP)
232d747510SLisandro Dalcin   #define PetscOptNameCmp(a, b) strcasecmp(a, b)
242d747510SLisandro Dalcin #elif defined(PETSC_HAVE_STRICMP)
252d747510SLisandro Dalcin   #define PetscOptNameCmp(a, b) stricmp(a, b)
262d747510SLisandro Dalcin #else
272d747510SLisandro Dalcin   #define PetscOptNameCmp(a, b) Error_strcasecmp_not_found
282d747510SLisandro Dalcin #endif
292d747510SLisandro Dalcin 
302d747510SLisandro Dalcin #include <petsc/private/hashtable.h>
312d747510SLisandro Dalcin 
322d747510SLisandro Dalcin /* This assumes ASCII encoding and ignores locale settings */
332d747510SLisandro Dalcin /* Using tolower() is about 2X slower in microbenchmarks   */
34d71ae5a4SJacob Faibussowitsch static inline int PetscToLower(int c)
35d71ae5a4SJacob Faibussowitsch {
362d747510SLisandro Dalcin   return ((c >= 'A') & (c <= 'Z')) ? c + 'a' - 'A' : c;
372d747510SLisandro Dalcin }
382d747510SLisandro Dalcin 
392d747510SLisandro Dalcin /* Bob Jenkins's one at a time hash function (case-insensitive) */
40d71ae5a4SJacob Faibussowitsch static inline unsigned int PetscOptHash(const char key[])
41d71ae5a4SJacob Faibussowitsch {
422d747510SLisandro Dalcin   unsigned int hash = 0;
432d747510SLisandro Dalcin   while (*key) {
442d747510SLisandro Dalcin     hash += PetscToLower(*key++);
452d747510SLisandro Dalcin     hash += hash << 10;
462d747510SLisandro Dalcin     hash ^= hash >> 6;
472d747510SLisandro Dalcin   }
482d747510SLisandro Dalcin   hash += hash << 3;
492d747510SLisandro Dalcin   hash ^= hash >> 11;
502d747510SLisandro Dalcin   hash += hash << 15;
512d747510SLisandro Dalcin   return hash;
522d747510SLisandro Dalcin }
532d747510SLisandro Dalcin 
54d71ae5a4SJacob Faibussowitsch static inline int PetscOptEqual(const char a[], const char b[])
55d71ae5a4SJacob Faibussowitsch {
562d747510SLisandro Dalcin   return !PetscOptNameCmp(a, b);
572d747510SLisandro Dalcin }
582d747510SLisandro Dalcin 
592d747510SLisandro Dalcin KHASH_INIT(HO, kh_cstr_t, int, 1, PetscOptHash, PetscOptEqual)
602d747510SLisandro Dalcin 
6174e0666dSJed Brown #define MAXPREFIXES        25
622d747510SLisandro Dalcin #define MAXOPTIONSMONITORS 5
63e5c89e4eSSatish Balay 
649355ec05SMatthew G. Knepley const char *PetscOptionSources[] = {"code", "command line", "file", "environment"};
659355ec05SMatthew G. Knepley 
669355ec05SMatthew G. Knepley // This table holds all the options set by the user
674416b707SBarry Smith struct _n_PetscOptions {
683de2bfdfSBarry Smith   PetscOptions previous;
699355ec05SMatthew G. Knepley 
702d747510SLisandro Dalcin   int                N;      /* number of options */
719355ec05SMatthew G. Knepley   int                Nalloc; /* number of allocated options */
729355ec05SMatthew G. Knepley   char             **names;  /* option names */
739355ec05SMatthew G. Knepley   char             **values; /* option values */
749355ec05SMatthew G. Knepley   PetscBool         *used;   /* flag option use */
759355ec05SMatthew G. Knepley   PetscOptionSource *source; /* source for option value */
76c5b5d8d5SVaclav Hapla   PetscBool          precedentProcessed;
77081c24baSBoyana Norris 
782d747510SLisandro Dalcin   /* Hash table */
792d747510SLisandro Dalcin   khash_t(HO) *ht;
802d747510SLisandro Dalcin 
812d747510SLisandro Dalcin   /* Prefixes */
822d747510SLisandro Dalcin   int  prefixind;
832d747510SLisandro Dalcin   int  prefixstack[MAXPREFIXES];
849355ec05SMatthew G. Knepley   char prefix[PETSC_MAX_OPTION_NAME];
852d747510SLisandro Dalcin 
862d747510SLisandro Dalcin   /* Aliases */
879355ec05SMatthew G. Knepley   int    Na;       /* number or aliases */
889355ec05SMatthew G. Knepley   int    Naalloc;  /* number of allocated aliases */
899355ec05SMatthew G. Knepley   char **aliases1; /* aliased */
909355ec05SMatthew G. Knepley   char **aliases2; /* aliasee */
912d747510SLisandro Dalcin 
922d747510SLisandro Dalcin   /* Help */
932d747510SLisandro Dalcin   PetscBool help;       /* flag whether "-help" is in the database */
94d314f959SVaclav Hapla   PetscBool help_intro; /* flag whether "-help intro" is in the database */
952d747510SLisandro Dalcin 
962d747510SLisandro Dalcin   /* Monitors */
97c5b5d8d5SVaclav Hapla   PetscBool monitorFromOptions, monitorCancel;
989355ec05SMatthew G. Knepley   PetscErrorCode (*monitor[MAXOPTIONSMONITORS])(const char[], const char[], PetscOptionSource, void *); /* returns control to user after */
999355ec05SMatthew G. Knepley   PetscErrorCode (*monitordestroy[MAXOPTIONSMONITORS])(void **);                                        /* callback for monitor destruction */
100081c24baSBoyana Norris   void    *monitorcontext[MAXOPTIONSMONITORS];                                                          /* to pass arbitrary user data into monitor */
101081c24baSBoyana Norris   PetscInt numbermonitors;                                                                              /* to, for instance, detect options being set */
1024416b707SBarry Smith };
103e5c89e4eSSatish Balay 
104b4205f0bSBarry Smith static PetscOptions defaultoptions = NULL; /* the options database routines query this object for options */
1052d747510SLisandro Dalcin 
106aaa8cc7dSPierre Jolivet /* list of options which precede others, i.e., are processed in PetscOptionsProcessPrecedentFlags() */
107660278c0SBarry Smith /* these options can only take boolean values, the code will crash if given a non-boolean value */
108660278c0SBarry Smith static const char *precedentOptions[] = {"-petsc_ci", "-options_monitor", "-options_monitor_cancel", "-help", "-skip_petscrc"};
1099371c9d4SSatish Balay enum PetscPrecedentOption {
1109371c9d4SSatish Balay   PO_CI_ENABLE,
1119371c9d4SSatish Balay   PO_OPTIONS_MONITOR,
1129371c9d4SSatish Balay   PO_OPTIONS_MONITOR_CANCEL,
1139371c9d4SSatish Balay   PO_HELP,
1149371c9d4SSatish Balay   PO_SKIP_PETSCRC,
1159371c9d4SSatish Balay   PO_NUM
1169371c9d4SSatish Balay };
117c5b5d8d5SVaclav Hapla 
1189355ec05SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscOptionsSetValue_Private(PetscOptions, const char[], const char[], int *, PetscOptionSource);
1199355ec05SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscOptionsInsertStringYAML_Private(PetscOptions, const char[], PetscOptionSource);
120c5b5d8d5SVaclav Hapla 
121081c24baSBoyana Norris /*
122081c24baSBoyana Norris     Options events monitor
123081c24baSBoyana Norris */
1249355ec05SMatthew G. Knepley static PetscErrorCode PetscOptionsMonitor(PetscOptions options, const char name[], const char value[], PetscOptionSource source)
125d71ae5a4SJacob Faibussowitsch {
126e5c89e4eSSatish Balay   PetscFunctionBegin;
127c5b5d8d5SVaclav Hapla   if (!value) value = "";
1289355ec05SMatthew G. Knepley   if (options->monitorFromOptions) PetscCall(PetscOptionsMonitorDefault(name, value, source, NULL));
1299355ec05SMatthew G. Knepley   for (PetscInt i = 0; i < options->numbermonitors; i++) PetscCall((*options->monitor[i])(name, value, source, options->monitorcontext[i]));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131e5c89e4eSSatish Balay }
132e5c89e4eSSatish Balay 
1332d747510SLisandro Dalcin /*@
1342d747510SLisandro Dalcin    PetscOptionsCreate - Creates an empty options database.
135e5c89e4eSSatish Balay 
13620f4b53cSBarry Smith    Logically Collective
1371c9f3c13SBarry Smith 
138e5c89e4eSSatish Balay    Output Parameter:
1392d747510SLisandro Dalcin .  options - Options database object
140e5c89e4eSSatish Balay 
141e5c89e4eSSatish Balay    Level: advanced
142e5c89e4eSSatish Balay 
143811af0c4SBarry Smith    Note:
144811af0c4SBarry Smith    Though PETSc has a concept of multiple options database the current code uses a single default `PetscOptions` object
145811af0c4SBarry Smith 
146811af0c4SBarry Smith    Developer Notes:
147811af0c4SBarry Smith    We may want eventually to pass a `MPI_Comm` to determine the ownership of the object
148811af0c4SBarry Smith 
149811af0c4SBarry Smith    This object never got developed after being introduced, it is not clear that supporting multiple `PetscOptions` objects is useful
1501c9f3c13SBarry Smith 
151db781477SPatrick Sanan .seealso: `PetscOptionsDestroy()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
152e5c89e4eSSatish Balay @*/
153d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsCreate(PetscOptions *options)
154d71ae5a4SJacob Faibussowitsch {
15539a651e2SJacob Faibussowitsch   PetscFunctionBegin;
15639a651e2SJacob Faibussowitsch   PetscValidPointer(options, 1);
1572d747510SLisandro Dalcin   *options = (PetscOptions)calloc(1, sizeof(**options));
15839a651e2SJacob Faibussowitsch   PetscCheck(*options, PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate the options database");
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1602d747510SLisandro Dalcin }
1612d747510SLisandro Dalcin 
1622d747510SLisandro Dalcin /*@
1632d747510SLisandro Dalcin     PetscOptionsDestroy - Destroys an option database.
1642d747510SLisandro Dalcin 
16520f4b53cSBarry Smith     Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
1661c9f3c13SBarry Smith 
1672d747510SLisandro Dalcin   Input Parameter:
168811af0c4SBarry Smith .  options - the `PetscOptions` object
1692d747510SLisandro Dalcin 
1703de2bfdfSBarry Smith    Level: advanced
1712d747510SLisandro Dalcin 
172db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
1732d747510SLisandro Dalcin @*/
174d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
175d71ae5a4SJacob Faibussowitsch {
176362febeeSStefano Zampini   PetscFunctionBegin;
177*10c654e6SJacob Faibussowitsch   PetscValidPointer(options, 1);
1783ba16761SJacob Faibussowitsch   if (!*options) PetscFunctionReturn(PETSC_SUCCESS);
1795f80ce2aSJacob Faibussowitsch   PetscCheck(!(*options)->previous, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "You are destroying an option that has been used with PetscOptionsPush() but does not have a corresponding PetscOptionsPop()");
1809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsClear(*options));
1812d747510SLisandro Dalcin   /* XXX what about monitors ? */
1822800570dSLisandro Dalcin   free(*options);
1832d747510SLisandro Dalcin   *options = NULL;
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
185e5c89e4eSSatish Balay }
186e5c89e4eSSatish Balay 
1872d747510SLisandro Dalcin /*
1882d747510SLisandro Dalcin     PetscOptionsCreateDefault - Creates the default global options database
1892d747510SLisandro Dalcin */
190d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsCreateDefault(void)
191d71ae5a4SJacob Faibussowitsch {
19239a651e2SJacob Faibussowitsch   PetscFunctionBegin;
1939566063dSJacob Faibussowitsch   if (PetscUnlikely(!defaultoptions)) PetscCall(PetscOptionsCreate(&defaultoptions));
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1952d747510SLisandro Dalcin }
1962d747510SLisandro Dalcin 
197b4205f0bSBarry Smith /*@
198811af0c4SBarry Smith       PetscOptionsPush - Push a new `PetscOptions` object as the default provider of options
1991c9f3c13SBarry Smith                          Allows using different parts of a code to use different options databases
200b4205f0bSBarry Smith 
201b4205f0bSBarry Smith   Logically Collective
202b4205f0bSBarry Smith 
203b4205f0bSBarry Smith   Input Parameter:
204811af0c4SBarry Smith .   opt - the options obtained with `PetscOptionsCreate()`
205b4205f0bSBarry Smith 
20620f4b53cSBarry Smith    Level: advanced
20720f4b53cSBarry Smith 
208b4205f0bSBarry Smith   Notes:
209811af0c4SBarry Smith   Use `PetscOptionsPop()` to return to the previous default options database
2101c9f3c13SBarry Smith 
211811af0c4SBarry Smith   The collectivity of this routine is complex; only the MPI ranks that call this routine will
2121c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
2131c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
2141c9f3c13SBarry Smith   on different ranks.
215b4205f0bSBarry Smith 
216811af0c4SBarry Smith   Developer Note:
217811af0c4SBarry Smith   Though this functionality has been provided it has never been used in PETSc and might be removed.
218811af0c4SBarry Smith 
219db781477SPatrick Sanan .seealso: `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
220b4205f0bSBarry Smith @*/
221d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPush(PetscOptions opt)
222d71ae5a4SJacob Faibussowitsch {
223b4205f0bSBarry Smith   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
225b4205f0bSBarry Smith   opt->previous  = defaultoptions;
226b4205f0bSBarry Smith   defaultoptions = opt;
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228b4205f0bSBarry Smith }
229b4205f0bSBarry Smith 
230b4205f0bSBarry Smith /*@
231811af0c4SBarry Smith       PetscOptionsPop - Pop the most recent `PetscOptionsPush()` to return to the previous default options
232b4205f0bSBarry Smith 
23320f4b53cSBarry Smith       Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
234b4205f0bSBarry Smith 
2353de2bfdfSBarry Smith    Level: advanced
2363de2bfdfSBarry Smith 
237db781477SPatrick Sanan .seealso: `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
238b4205f0bSBarry Smith @*/
239d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPop(void)
240d71ae5a4SJacob Faibussowitsch {
2413de2bfdfSBarry Smith   PetscOptions current = defaultoptions;
2423de2bfdfSBarry Smith 
243b4205f0bSBarry Smith   PetscFunctionBegin;
24428b400f6SJacob Faibussowitsch   PetscCheck(defaultoptions, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing default options");
24528b400f6SJacob Faibussowitsch   PetscCheck(defaultoptions->previous, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscOptionsPop() called too many times");
246b4205f0bSBarry Smith   defaultoptions    = defaultoptions->previous;
2473de2bfdfSBarry Smith   current->previous = NULL;
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
249b4205f0bSBarry Smith }
250b4205f0bSBarry Smith 
2512d747510SLisandro Dalcin /*
2522d747510SLisandro Dalcin     PetscOptionsDestroyDefault - Destroys the default global options database
2532d747510SLisandro Dalcin */
254d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDestroyDefault(void)
255d71ae5a4SJacob Faibussowitsch {
25639a651e2SJacob Faibussowitsch   PetscFunctionBegin;
2573ba16761SJacob Faibussowitsch   if (!defaultoptions) PetscFunctionReturn(PETSC_SUCCESS);
2583de2bfdfSBarry Smith   /* Destroy any options that the user forgot to pop */
2593de2bfdfSBarry Smith   while (defaultoptions->previous) {
26039a651e2SJacob Faibussowitsch     PetscOptions tmp = defaultoptions;
26139a651e2SJacob Faibussowitsch 
2629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPop());
2639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsDestroy(&tmp));
2643de2bfdfSBarry Smith   }
2659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroy(&defaultoptions));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267e5c89e4eSSatish Balay }
268e5c89e4eSSatish Balay 
26994ef8ddeSSatish Balay /*@C
2707cd08cecSJed Brown    PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.
2713fc1eb6aSBarry Smith 
27220f4b53cSBarry Smith    Not Collective
2731c9f3c13SBarry Smith 
2743fc1eb6aSBarry Smith    Input Parameter:
2752d747510SLisandro Dalcin .  key - string to check if valid
2763fc1eb6aSBarry Smith 
2773fc1eb6aSBarry Smith    Output Parameter:
278811af0c4SBarry Smith .  valid - `PETSC_TRUE` if a valid key
2793fc1eb6aSBarry Smith 
280f6680f47SSatish Balay    Level: intermediate
2813fc1eb6aSBarry Smith @*/
282d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsValidKey(const char key[], PetscBool *valid)
283d71ae5a4SJacob Faibussowitsch {
284f603b5e9SToby Isaac   char *ptr;
2857c5db45bSBarry Smith 
28696fc60bcSBarry Smith   PetscFunctionBegin;
2872d747510SLisandro Dalcin   if (key) PetscValidCharPointer(key, 1);
288dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(valid, 2);
2892d747510SLisandro Dalcin   *valid = PETSC_FALSE;
2903ba16761SJacob Faibussowitsch   if (!key) PetscFunctionReturn(PETSC_SUCCESS);
2913ba16761SJacob Faibussowitsch   if (key[0] != '-') PetscFunctionReturn(PETSC_SUCCESS);
2922d747510SLisandro Dalcin   if (key[1] == '-') key++;
2933ba16761SJacob Faibussowitsch   if (!isalpha((int)key[1])) PetscFunctionReturn(PETSC_SUCCESS);
2942d747510SLisandro Dalcin   (void)strtod(key, &ptr);
2953ba16761SJacob Faibussowitsch   if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) PetscFunctionReturn(PETSC_SUCCESS);
2962d747510SLisandro Dalcin   *valid = PETSC_TRUE;
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29896fc60bcSBarry Smith }
29996fc60bcSBarry Smith 
300*10c654e6SJacob Faibussowitsch static PetscErrorCode PetscOptionsInsertString_Private(PetscOptions options, const char in_str[], PetscOptionSource source)
301d71ae5a4SJacob Faibussowitsch {
302d06005cbSLisandro Dalcin   char      *first, *second;
3039c9d3cfdSBarry Smith   PetscToken token;
304e5c89e4eSSatish Balay 
305e5c89e4eSSatish Balay   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(in_str, ' ', &token));
3079566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &first));
30896fc60bcSBarry Smith   while (first) {
309d06005cbSLisandro Dalcin     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
310*10c654e6SJacob Faibussowitsch 
3119566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_file", &isfile));
3129566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_file_yaml", &isfileyaml));
3139566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_string_yaml", &isstringyaml));
3149566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-prefix_push", &ispush));
3159566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-prefix_pop", &ispop));
3169566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(first, &key));
317d06005cbSLisandro Dalcin     if (!key) {
3189566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
319d06005cbSLisandro Dalcin     } else if (isfile) {
3209566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
321*10c654e6SJacob Faibussowitsch       PetscCall(PetscOptionsInsertFile(PETSC_COMM_SELF, options, second, PETSC_TRUE));
3229566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
323d06005cbSLisandro Dalcin     } else if (isfileyaml) {
3249566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
325*10c654e6SJacob Faibussowitsch       PetscCall(PetscOptionsInsertFileYAML(PETSC_COMM_SELF, options, second, PETSC_TRUE));
3269566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
327d06005cbSLisandro Dalcin     } else if (isstringyaml) {
3289566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3299355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, second, source));
3309566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
331d06005cbSLisandro Dalcin     } else if (ispush) {
3329566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3339566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPush(options, second));
3349566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
3359db968c8SJed Brown     } else if (ispop) {
3369566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPop(options));
3379566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
338d06005cbSLisandro Dalcin     } else {
3399566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3409566063dSJacob Faibussowitsch       PetscCall(PetscOptionsValidKey(second, &key));
34196fc60bcSBarry Smith       if (!key) {
3429355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, first, second, NULL, source));
3439566063dSJacob Faibussowitsch         PetscCall(PetscTokenFind(token, &first));
34496fc60bcSBarry Smith       } else {
3459355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, first, NULL, NULL, source));
34696fc60bcSBarry Smith         first = second;
34796fc60bcSBarry Smith       }
348e5c89e4eSSatish Balay     }
349e5c89e4eSSatish Balay   }
3509566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
352e5c89e4eSSatish Balay }
353e5c89e4eSSatish Balay 
3549355ec05SMatthew G. Knepley /*@C
3559355ec05SMatthew G. Knepley    PetscOptionsInsertString - Inserts options into the database from a string
3569355ec05SMatthew G. Knepley 
3579355ec05SMatthew G. Knepley    Logically Collective
3589355ec05SMatthew G. Knepley 
3599355ec05SMatthew G. Knepley    Input Parameters:
3609355ec05SMatthew G. Knepley +  options - options object
3619355ec05SMatthew G. Knepley -  in_str - string that contains options separated by blanks
3629355ec05SMatthew G. Knepley 
3639355ec05SMatthew G. Knepley    Level: intermediate
3649355ec05SMatthew G. Knepley 
36520f4b53cSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
3669355ec05SMatthew G. Knepley   have the affect of these options. If some processes that create objects call this routine and others do
3679355ec05SMatthew G. Knepley   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
3689355ec05SMatthew G. Knepley   on different ranks.
3699355ec05SMatthew G. Knepley 
3709355ec05SMatthew G. Knepley    Contributed by Boyana Norris
3719355ec05SMatthew G. Knepley 
3729355ec05SMatthew G. Knepley .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
3739355ec05SMatthew G. Knepley           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3749355ec05SMatthew G. Knepley           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3759355ec05SMatthew G. Knepley           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3769355ec05SMatthew G. Knepley           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3779355ec05SMatthew G. Knepley           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsInsertFile()`
3789355ec05SMatthew G. Knepley @*/
3799355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsInsertString(PetscOptions options, const char in_str[])
3809355ec05SMatthew G. Knepley {
3819355ec05SMatthew G. Knepley   PetscFunctionBegin;
3829355ec05SMatthew G. Knepley   PetscCall(PetscOptionsInsertString_Private(options, in_str, PETSC_OPT_CODE));
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3849355ec05SMatthew G. Knepley }
3859355ec05SMatthew G. Knepley 
3863fc1eb6aSBarry Smith /*
3873fc1eb6aSBarry Smith     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
3883fc1eb6aSBarry Smith */
389d71ae5a4SJacob Faibussowitsch static char *Petscgetline(FILE *f)
390d71ae5a4SJacob Faibussowitsch {
3915fa91da5SBarry Smith   size_t size = 0;
3925fa91da5SBarry Smith   size_t len  = 0;
3935fa91da5SBarry Smith   size_t last = 0;
3940298fd71SBarry Smith   char  *buf  = NULL;
3955fa91da5SBarry Smith 
39602c9f0b5SLisandro Dalcin   if (feof(f)) return NULL;
3975fa91da5SBarry Smith   do {
3985fa91da5SBarry Smith     size += 1024;                             /* BUFSIZ is defined as "the optimal read size for this platform" */
3996e0c8459SSatish Balay     buf = (char *)realloc((void *)buf, size); /* realloc(NULL,n) is the same as malloc(n) */
4005fa91da5SBarry Smith     /* Actually do the read. Note that fgets puts a terminal '\0' on the
4015fa91da5SBarry Smith     end of the string, so we make sure we overwrite this */
402e86f3e45SDave May     if (!fgets(buf + len, 1024, f)) buf[len] = 0;
4033ba16761SJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, PetscStrlen(buf, &len));
4045fa91da5SBarry Smith     last = len - 1;
4055fa91da5SBarry Smith   } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
40608ac41f7SSatish Balay   if (len) return buf;
4075fa91da5SBarry Smith   free(buf);
40802c9f0b5SLisandro Dalcin   return NULL;
4095fa91da5SBarry Smith }
4105fa91da5SBarry Smith 
411d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsFilename(MPI_Comm comm, const char file[], char filename[PETSC_MAX_PATH_LEN], PetscBool *yaml)
412d71ae5a4SJacob Faibussowitsch {
413be10d61cSLisandro Dalcin   char fname[PETSC_MAX_PATH_LEN + 8], path[PETSC_MAX_PATH_LEN + 8], *tail;
414e5c89e4eSSatish Balay 
415be10d61cSLisandro Dalcin   PetscFunctionBegin;
416362febeeSStefano Zampini   *yaml = PETSC_FALSE;
4179566063dSJacob Faibussowitsch   PetscCall(PetscStrreplace(comm, file, fname, sizeof(fname)));
4189566063dSJacob Faibussowitsch   PetscCall(PetscFixFilename(fname, path));
4199566063dSJacob Faibussowitsch   PetscCall(PetscStrendswith(path, ":yaml", yaml));
420be10d61cSLisandro Dalcin   if (*yaml) {
4219566063dSJacob Faibussowitsch     PetscCall(PetscStrrchr(path, ':', &tail));
422be10d61cSLisandro Dalcin     tail[-1] = 0; /* remove ":yaml" suffix from path */
423be10d61cSLisandro Dalcin   }
4249566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename, path, PETSC_MAX_PATH_LEN));
425a1d2f846SLisandro Dalcin   /* check for standard YAML and JSON filename extensions */
4269566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yaml", yaml));
4279566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yml", yaml));
4289566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".json", yaml));
429a1d2f846SLisandro Dalcin   if (!*yaml) { /* check file contents */
430a1d2f846SLisandro Dalcin     PetscMPIInt rank;
4319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
432dd400576SPatrick Sanan     if (rank == 0) {
433a1d2f846SLisandro Dalcin       FILE *fh = fopen(filename, "r");
434a1d2f846SLisandro Dalcin       if (fh) {
435a1d2f846SLisandro Dalcin         char buf[6] = "";
436a1d2f846SLisandro Dalcin         if (fread(buf, 1, 6, fh) > 0) {
4379566063dSJacob Faibussowitsch           PetscCall(PetscStrncmp(buf, "%YAML ", 6, yaml));          /* check for '%YAML' tag */
4389566063dSJacob Faibussowitsch           if (!*yaml) PetscCall(PetscStrncmp(buf, "---", 3, yaml)); /* check for document start */
439a1d2f846SLisandro Dalcin         }
440a1d2f846SLisandro Dalcin         (void)fclose(fh);
441a1d2f846SLisandro Dalcin       }
442a1d2f846SLisandro Dalcin     }
4439566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(yaml, 1, MPIU_BOOL, 0, comm));
444a1d2f846SLisandro Dalcin   }
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446be10d61cSLisandro Dalcin }
447e5c89e4eSSatish Balay 
448d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsInsertFilePetsc(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
449d71ae5a4SJacob Faibussowitsch {
4508c0b561eSLisandro Dalcin   char       *string, *vstring = NULL, *astring = NULL, *packed = NULL;
4517fb43599SVaclav Hapla   char       *tokens[4];
45213e3f751SJed Brown   size_t      i, len, bytes;
453e5c89e4eSSatish Balay   FILE       *fd;
4547fb43599SVaclav Hapla   PetscToken  token = NULL;
455ed9cf6e9SBarry Smith   int         err;
456bbcf679cSJacob Faibussowitsch   char       *cmatch = NULL;
457581bbe83SVaclav Hapla   const char  cmt    = '#';
4589210b8eaSVaclav Hapla   PetscInt    line   = 1;
4593a018368SJed Brown   PetscMPIInt rank, cnt = 0, acnt = 0, counts[2];
4609210b8eaSVaclav Hapla   PetscBool   isdir, alias = PETSC_FALSE, valid;
461e5c89e4eSSatish Balay 
462e5c89e4eSSatish Balay   PetscFunctionBegin;
4639566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(tokens, sizeof(tokens)));
4649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
465dd400576SPatrick Sanan   if (rank == 0) {
4668c0b561eSLisandro Dalcin     char fpath[PETSC_MAX_PATH_LEN];
4678c0b561eSLisandro Dalcin     char fname[PETSC_MAX_PATH_LEN];
46805c7dedfSBarry Smith 
4699566063dSJacob Faibussowitsch     PetscCall(PetscStrreplace(PETSC_COMM_SELF, file, fpath, sizeof(fpath)));
4709566063dSJacob Faibussowitsch     PetscCall(PetscFixFilename(fpath, fname));
4718c0b561eSLisandro Dalcin 
472e5c89e4eSSatish Balay     fd = fopen(fname, "r");
4739566063dSJacob Faibussowitsch     PetscCall(PetscTestDirectory(fname, 'r', &isdir));
47408401ef6SPierre Jolivet     PetscCheck(!isdir || !require, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified options file %s is a directory", fname);
475ad38b122SPatrick Sanan     if (fd && !isdir) {
4763a018368SJed Brown       PetscSegBuffer vseg, aseg;
4779566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferCreate(1, 4000, &vseg));
4789566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferCreate(1, 2000, &aseg));
4793a018368SJed Brown 
4809b754dc9SBarry Smith       /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
4819566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Opened options file %s\n", file));
482e24ecc5dSJed Brown 
4835fa91da5SBarry Smith       while ((string = Petscgetline(fd))) {
4844704e885SBarry Smith         /* eliminate comments from each line */
4859566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(string, cmt, &cmatch));
48690f79514SSatish Balay         if (cmatch) *cmatch = 0;
4879566063dSJacob Faibussowitsch         PetscCall(PetscStrlen(string, &len));
4885981331cSSatish Balay         /* replace tabs, ^M, \n with " " */
489e5c89e4eSSatish Balay         for (i = 0; i < len; i++) {
490ad540459SPierre Jolivet           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') string[i] = ' ';
491e5c89e4eSSatish Balay         }
4929566063dSJacob Faibussowitsch         PetscCall(PetscTokenCreate(string, ' ', &token));
4939566063dSJacob Faibussowitsch         PetscCall(PetscTokenFind(token, &tokens[0]));
4947fb43599SVaclav Hapla         if (!tokens[0]) {
49502b0d46eSSatish Balay           goto destroy;
4967fb43599SVaclav Hapla         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
4979566063dSJacob Faibussowitsch           PetscCall(PetscTokenFind(token, &tokens[0]));
49890f79514SSatish Balay         }
49948a46eb9SPierre Jolivet         for (i = 1; i < 4; i++) PetscCall(PetscTokenFind(token, &tokens[i]));
5007fb43599SVaclav Hapla         if (!tokens[0]) {
5012662f744SSatish Balay           goto destroy;
5027fb43599SVaclav Hapla         } else if (tokens[0][0] == '-') {
5039566063dSJacob Faibussowitsch           PetscCall(PetscOptionsValidKey(tokens[0], &valid));
50428b400f6SJacob Faibussowitsch           PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": invalid option %s", fname, line, tokens[0]);
5059566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(tokens[0], &len));
5069566063dSJacob Faibussowitsch           PetscCall(PetscSegBufferGet(vseg, len + 1, &vstring));
5079566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(vstring, tokens[0], len));
508e24ecc5dSJed Brown           vstring[len] = ' ';
5097fb43599SVaclav Hapla           if (tokens[1]) {
5109566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
51128b400f6SJacob Faibussowitsch             PetscCheck(!valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": cannot specify two options per line (%s %s)", fname, line, tokens[0], tokens[1]);
5129566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[1], &len));
5139566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(vseg, len + 3, &vstring));
514e24ecc5dSJed Brown             vstring[0] = '"';
5159566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(vstring + 1, tokens[1], len));
516e24ecc5dSJed Brown             vstring[len + 1] = '"';
517e24ecc5dSJed Brown             vstring[len + 2] = ' ';
51809192fe3SBarry Smith           }
51990f79514SSatish Balay         } else {
5209566063dSJacob Faibussowitsch           PetscCall(PetscStrcasecmp(tokens[0], "alias", &alias));
5219210b8eaSVaclav Hapla           if (alias) {
5229566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
52328b400f6SJacob Faibussowitsch             PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": invalid aliased option %s", fname, line, tokens[1]);
52408401ef6SPierre Jolivet             PetscCheck(tokens[2], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": alias missing for %s", fname, line, tokens[1]);
5259566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[2], &valid));
52628b400f6SJacob Faibussowitsch             PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": invalid aliasee option %s", fname, line, tokens[2]);
5279566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[1], &len));
5289566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
5299566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(astring, tokens[1], len));
530e24ecc5dSJed Brown             astring[len] = ' ';
531e24ecc5dSJed Brown 
5329566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[2], &len));
5339566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
5349566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(astring, tokens[2], len));
535e24ecc5dSJed Brown             astring[len] = ' ';
53698921bdaSJacob Faibussowitsch           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown first token in options file %s line %" PetscInt_FMT ": %s", fname, line, tokens[0]);
5379210b8eaSVaclav Hapla         }
5389210b8eaSVaclav Hapla         {
5399210b8eaSVaclav Hapla           const char *extraToken = alias ? tokens[3] : tokens[2];
54028b400f6SJacob Faibussowitsch           PetscCheck(!extraToken, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": extra token %s", fname, line, extraToken);
541e5c89e4eSSatish Balay         }
54202b0d46eSSatish Balay       destroy:
5434b40f50bSBarry Smith         free(string);
5449566063dSJacob Faibussowitsch         PetscCall(PetscTokenDestroy(&token));
5459210b8eaSVaclav Hapla         alias = PETSC_FALSE;
5469210b8eaSVaclav Hapla         line++;
547e5c89e4eSSatish Balay       }
548ed9cf6e9SBarry Smith       err = fclose(fd);
54928b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file %s", fname);
5509566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGetSize(aseg, &bytes)); /* size without null termination */
5519566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(bytes, &acnt));
5529566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(aseg, 1, &astring));
553e24ecc5dSJed Brown       astring[0] = 0;
5549566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGetSize(vseg, &bytes)); /* size without null termination */
5559566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(bytes, &cnt));
5569566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(vseg, 1, &vstring));
557e24ecc5dSJed Brown       vstring[0] = 0;
5589566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
5599566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferExtractTo(aseg, packed));
5609566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferExtractTo(vseg, packed + acnt + 1));
5619566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferDestroy(&aseg));
5629566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferDestroy(&vseg));
56328b400f6SJacob Faibussowitsch     } else PetscCheck(!require, PETSC_COMM_SELF, PETSC_ERR_USER, "Unable to open options file %s", fname);
5649b754dc9SBarry Smith   }
56505c7dedfSBarry Smith 
5663a018368SJed Brown   counts[0] = acnt;
5673a018368SJed Brown   counts[1] = cnt;
5684201f521SBarry Smith   err       = MPI_Bcast(counts, 2, MPI_INT, 0, comm);
56928b400f6SJacob Faibussowitsch   PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in first MPI collective call, could be caused by using an incorrect mpiexec or a network problem, it can be caused by having VPN running: see https://petsc.org/release/faq/");
5703a018368SJed Brown   acnt = counts[0];
5713a018368SJed Brown   cnt  = counts[1];
57248a46eb9SPierre Jolivet   if (rank) PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
5733a018368SJed Brown   if (acnt || cnt) {
5749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(packed, 2 + acnt + cnt, MPI_CHAR, 0, comm));
5753a018368SJed Brown     astring = packed;
5763a018368SJed Brown     vstring = packed + acnt + 1;
5773a018368SJed Brown   }
5783a018368SJed Brown 
5799b754dc9SBarry Smith   if (acnt) {
5809566063dSJacob Faibussowitsch     PetscCall(PetscTokenCreate(astring, ' ', &token));
5819566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &tokens[0]));
5827fb43599SVaclav Hapla     while (tokens[0]) {
5839566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &tokens[1]));
5849566063dSJacob Faibussowitsch       PetscCall(PetscOptionsSetAlias(options, tokens[0], tokens[1]));
5859566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &tokens[0]));
5869b754dc9SBarry Smith     }
5879566063dSJacob Faibussowitsch     PetscCall(PetscTokenDestroy(&token));
5889b754dc9SBarry Smith   }
5899b754dc9SBarry Smith 
5909355ec05SMatthew G. Knepley   if (cnt) PetscCall(PetscOptionsInsertString_Private(options, vstring, PETSC_OPT_FILE));
5919566063dSJacob Faibussowitsch   PetscCall(PetscFree(packed));
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593e5c89e4eSSatish Balay }
594e5c89e4eSSatish Balay 
595d06005cbSLisandro Dalcin /*@C
596be10d61cSLisandro Dalcin      PetscOptionsInsertFile - Inserts options into the database from a file.
597be10d61cSLisandro Dalcin 
598be10d61cSLisandro Dalcin      Collective
599be10d61cSLisandro Dalcin 
600d8d19677SJose E. Roman   Input Parameters:
601811af0c4SBarry Smith +   comm - the processes that will share the options (usually `PETSC_COMM_WORLD`)
60220f4b53cSBarry Smith .   options - options database, use `NULL` for default global database
603be10d61cSLisandro Dalcin .   file - name of file,
604be10d61cSLisandro Dalcin            ".yml" and ".yaml" filename extensions are inserted as YAML options,
605be10d61cSLisandro Dalcin            append ":yaml" to filename to force YAML options.
606811af0c4SBarry Smith -   require - if `PETSC_TRUE` will generate an error if the file does not exist
607be10d61cSLisandro Dalcin 
60820f4b53cSBarry Smith   Level: developer
60920f4b53cSBarry Smith 
610be10d61cSLisandro Dalcin   Notes:
611be10d61cSLisandro Dalcin    Use  # for lines that are comments and which should be ignored.
612811af0c4SBarry Smith    Usually, instead of using this command, one should list the file name in the call to `PetscInitialize()`, this insures that certain options
61321532e8aSBarry Smith    such as `-log_view` or `-malloc_debug` are processed properly. This routine only sets options into the options database that will be processed by later
61421532e8aSBarry Smith    calls to `XXXSetFromOptions()`, it should not be used for options listed under PetscInitialize().
61521532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes in comm will
61621532e8aSBarry Smith    have the effect of these options. If some processes that create objects call this routine and others do
617be10d61cSLisandro Dalcin    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
618be10d61cSLisandro Dalcin    on different ranks.
619be10d61cSLisandro Dalcin 
620db781477SPatrick Sanan .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
621db781477SPatrick Sanan           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
622db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
623c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
624db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
625db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
626be10d61cSLisandro Dalcin @*/
627d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
628d71ae5a4SJacob Faibussowitsch {
629be10d61cSLisandro Dalcin   char      filename[PETSC_MAX_PATH_LEN];
630be10d61cSLisandro Dalcin   PetscBool yaml;
631be10d61cSLisandro Dalcin 
632be10d61cSLisandro Dalcin   PetscFunctionBegin;
6339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFilename(comm, file, filename, &yaml));
634be10d61cSLisandro Dalcin   if (yaml) {
6359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFileYAML(comm, options, filename, require));
636be10d61cSLisandro Dalcin   } else {
6379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFilePetsc(comm, options, filename, require));
638be10d61cSLisandro Dalcin   }
6393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
640be10d61cSLisandro Dalcin }
641be10d61cSLisandro Dalcin 
642be10d61cSLisandro Dalcin /*@C
643d06005cbSLisandro Dalcin    PetscOptionsInsertArgs - Inserts options into the database from a array of strings
644d06005cbSLisandro Dalcin 
645d06005cbSLisandro Dalcin    Logically Collective
646d06005cbSLisandro Dalcin 
647d8d19677SJose E. Roman    Input Parameters:
648d06005cbSLisandro Dalcin +  options - options object
6496aad120cSJose E. Roman .  argc - the array length
650d06005cbSLisandro Dalcin -  args - the string array
651d06005cbSLisandro Dalcin 
652d06005cbSLisandro Dalcin    Level: intermediate
653d06005cbSLisandro Dalcin 
654db781477SPatrick Sanan .seealso: `PetscOptions`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`
655d06005cbSLisandro Dalcin @*/
656d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsertArgs(PetscOptions options, int argc, char *args[])
657d71ae5a4SJacob Faibussowitsch {
658d06005cbSLisandro Dalcin   MPI_Comm     comm  = PETSC_COMM_WORLD;
659d06005cbSLisandro Dalcin   int          left  = PetscMax(argc, 0);
660d06005cbSLisandro Dalcin   char *const *eargs = args;
66185079163SJed Brown 
66285079163SJed Brown   PetscFunctionBegin;
66385079163SJed Brown   while (left) {
664d06005cbSLisandro Dalcin     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
6659566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_file", &isfile));
6669566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_file_yaml", &isfileyaml));
6679566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_string_yaml", &isstringyaml));
6689566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_push", &ispush));
6699566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_pop", &ispop));
6709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(eargs[0], &key));
671093de6efSBarry Smith     if (!key) {
6729371c9d4SSatish Balay       eargs++;
6739371c9d4SSatish Balay       left--;
674d06005cbSLisandro Dalcin     } else if (isfile) {
675cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file filename option");
6769566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInsertFile(comm, options, eargs[1], PETSC_TRUE));
6779371c9d4SSatish Balay       eargs += 2;
6789371c9d4SSatish Balay       left -= 2;
679d06005cbSLisandro Dalcin     } else if (isfileyaml) {
680cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file_yaml filename option");
6819566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInsertFileYAML(comm, options, eargs[1], PETSC_TRUE));
6829371c9d4SSatish Balay       eargs += 2;
6839371c9d4SSatish Balay       left -= 2;
684d06005cbSLisandro Dalcin     } else if (isstringyaml) {
685cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing string for -options_string_yaml string option");
6869355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, eargs[1], PETSC_OPT_CODE));
6879371c9d4SSatish Balay       eargs += 2;
6889371c9d4SSatish Balay       left -= 2;
689d06005cbSLisandro Dalcin     } else if (ispush) {
69008401ef6SPierre Jolivet       PetscCheck(left > 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option");
691cc73adaaSBarry Smith       PetscCheck(eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option (prefixes cannot start with '-')");
6929566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPush(options, eargs[1]));
6939371c9d4SSatish Balay       eargs += 2;
6949371c9d4SSatish Balay       left -= 2;
695d06005cbSLisandro Dalcin     } else if (ispop) {
6969566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPop(options));
6979371c9d4SSatish Balay       eargs++;
6989371c9d4SSatish Balay       left--;
6997935c3d8SJed Brown     } else {
7007935c3d8SJed Brown       PetscBool nextiskey = PETSC_FALSE;
7019566063dSJacob Faibussowitsch       if (left >= 2) PetscCall(PetscOptionsValidKey(eargs[1], &nextiskey));
70298b6bf53SJed Brown       if (left < 2 || nextiskey) {
7039355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], NULL, NULL, PETSC_OPT_COMMAND_LINE));
7049371c9d4SSatish Balay         eargs++;
7059371c9d4SSatish Balay         left--;
70685079163SJed Brown       } else {
7079355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], eargs[1], NULL, PETSC_OPT_COMMAND_LINE));
7089371c9d4SSatish Balay         eargs += 2;
7099371c9d4SSatish Balay         left -= 2;
71085079163SJed Brown       }
71185079163SJed Brown     }
7127935c3d8SJed Brown   }
7133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71485079163SJed Brown }
71585079163SJed Brown 
716*10c654e6SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt, const char *val[], const PetscBool set[], PetscBool *flg)
717d71ae5a4SJacob Faibussowitsch {
718c5b5d8d5SVaclav Hapla   PetscFunctionBegin;
719c5b5d8d5SVaclav Hapla   if (set[opt]) {
7209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(val[opt], flg));
721c5b5d8d5SVaclav Hapla   } else *flg = PETSC_FALSE;
7223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
723c5b5d8d5SVaclav Hapla }
724c5b5d8d5SVaclav Hapla 
725660278c0SBarry Smith /* Process options with absolute precedence, these are only processed from the command line, not the environment or files */
726d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options, int argc, char *args[], PetscBool *skip_petscrc, PetscBool *skip_petscrc_set)
727d71ae5a4SJacob Faibussowitsch {
728c5b5d8d5SVaclav Hapla   const char *const *opt = precedentOptions;
729c5b5d8d5SVaclav Hapla   const size_t       n   = PO_NUM;
730c5b5d8d5SVaclav Hapla   size_t             o;
731c5b5d8d5SVaclav Hapla   int                a;
732c5b5d8d5SVaclav Hapla   const char       **val;
7330c99d500SBarry Smith   char             **cval;
734660278c0SBarry Smith   PetscBool         *set, unneeded;
735c5b5d8d5SVaclav Hapla 
736c5b5d8d5SVaclav Hapla   PetscFunctionBegin;
7370c99d500SBarry Smith   PetscCall(PetscCalloc2(n, &cval, n, &set));
7380c99d500SBarry Smith   val = (const char **)cval;
739c5b5d8d5SVaclav Hapla 
740c5b5d8d5SVaclav Hapla   /* Look for options possibly set using PetscOptionsSetValue beforehand */
74148a46eb9SPierre Jolivet   for (o = 0; o < n; o++) PetscCall(PetscOptionsFindPair(options, NULL, opt[o], &val[o], &set[o]));
742c5b5d8d5SVaclav Hapla 
743a5b23f4aSJose E. Roman   /* Loop through all args to collect last occurring value of each option */
744c5b5d8d5SVaclav Hapla   for (a = 1; a < argc; a++) {
745c5b5d8d5SVaclav Hapla     PetscBool valid, eq;
746c5b5d8d5SVaclav Hapla 
7479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(args[a], &valid));
748c5b5d8d5SVaclav Hapla     if (!valid) continue;
749c5b5d8d5SVaclav Hapla     for (o = 0; o < n; o++) {
7509566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(args[a], opt[o], &eq));
751c5b5d8d5SVaclav Hapla       if (eq) {
752c5b5d8d5SVaclav Hapla         set[o] = PETSC_TRUE;
753c5b5d8d5SVaclav Hapla         if (a == argc - 1 || !args[a + 1] || !args[a + 1][0] || args[a + 1][0] == '-') val[o] = NULL;
754c5b5d8d5SVaclav Hapla         else val[o] = args[a + 1];
755c5b5d8d5SVaclav Hapla         break;
756c5b5d8d5SVaclav Hapla       }
757c5b5d8d5SVaclav Hapla     }
758c5b5d8d5SVaclav Hapla   }
759c5b5d8d5SVaclav Hapla 
760c5b5d8d5SVaclav Hapla   /* Process flags */
7619566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro));
762d314f959SVaclav Hapla   if (options->help_intro) options->help = PETSC_TRUE;
7639566063dSJacob Faibussowitsch   else PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_HELP, val, set, &options->help));
764660278c0SBarry Smith   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_CI_ENABLE, val, set, &unneeded));
765660278c0SBarry Smith   /* need to manage PO_CI_ENABLE option before the PetscOptionsMonitor is turned on, so its setting is not monitored */
7669355ec05SMatthew G. Knepley   if (set[PO_CI_ENABLE]) PetscCall(PetscOptionsSetValue_Private(options, opt[PO_CI_ENABLE], val[PO_CI_ENABLE], &a, PETSC_OPT_COMMAND_LINE));
7679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL, val, set, &options->monitorCancel));
7689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR, val, set, &options->monitorFromOptions));
7699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC, val, set, skip_petscrc));
770c5b5d8d5SVaclav Hapla   *skip_petscrc_set = set[PO_SKIP_PETSCRC];
771c5b5d8d5SVaclav Hapla 
772c5b5d8d5SVaclav Hapla   /* Store precedent options in database and mark them as used */
773660278c0SBarry Smith   for (o = 1; o < n; o++) {
774c5b5d8d5SVaclav Hapla     if (set[o]) {
7759355ec05SMatthew G. Knepley       PetscCall(PetscOptionsSetValue_Private(options, opt[o], val[o], &a, PETSC_OPT_COMMAND_LINE));
776d06005cbSLisandro Dalcin       options->used[a] = PETSC_TRUE;
777c5b5d8d5SVaclav Hapla     }
778c5b5d8d5SVaclav Hapla   }
7790c99d500SBarry Smith   PetscCall(PetscFree2(cval, set));
780c5b5d8d5SVaclav Hapla   options->precedentProcessed = PETSC_TRUE;
7813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
782c5b5d8d5SVaclav Hapla }
783c5b5d8d5SVaclav Hapla 
784d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options, const char name[], PetscBool *flg)
785d71ae5a4SJacob Faibussowitsch {
78639a651e2SJacob Faibussowitsch   PetscFunctionBegin;
78739a651e2SJacob Faibussowitsch   PetscValidBoolPointer(flg, 3);
788c5b5d8d5SVaclav Hapla   *flg = PETSC_FALSE;
789c5b5d8d5SVaclav Hapla   if (options->precedentProcessed) {
79039a651e2SJacob Faibussowitsch     for (int i = 0; i < PO_NUM; ++i) {
791c5b5d8d5SVaclav Hapla       if (!PetscOptNameCmp(precedentOptions[i], name)) {
792c5b5d8d5SVaclav Hapla         /* check if precedent option has been set already */
7939566063dSJacob Faibussowitsch         PetscCall(PetscOptionsFindPair(options, NULL, name, NULL, flg));
794c5b5d8d5SVaclav Hapla         if (*flg) break;
795c5b5d8d5SVaclav Hapla       }
796c5b5d8d5SVaclav Hapla     }
797c5b5d8d5SVaclav Hapla   }
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
799c5b5d8d5SVaclav Hapla }
80085079163SJed Brown 
801e5c89e4eSSatish Balay /*@C
802e5c89e4eSSatish Balay    PetscOptionsInsert - Inserts into the options database from the command line,
803e5c89e4eSSatish Balay                         the environmental variable and a file.
804e5c89e4eSSatish Balay 
805811af0c4SBarry Smith    Collective on `PETSC_COMM_WORLD`
8061c9f3c13SBarry Smith 
807e5c89e4eSSatish Balay    Input Parameters:
80820f4b53cSBarry Smith +  options - options database or `NULL` for the default global database
809c5929fdfSBarry Smith .  argc - count of number of command line arguments
810e5c89e4eSSatish Balay .  args - the command line arguments
811be10d61cSLisandro Dalcin -  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
81220f4b53cSBarry Smith           Use `NULL` or empty string to not check for code specific file.
813be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
814c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
815e5c89e4eSSatish Balay 
816081c24baSBoyana Norris    Options Database Keys:
817d06005cbSLisandro Dalcin +   -options_file <filename> - read options from a file
818d06005cbSLisandro Dalcin -   -options_file_yaml <filename> - read options from a YAML file
819c5b5d8d5SVaclav Hapla 
82020f4b53cSBarry Smith    Level: advanced
82120f4b53cSBarry Smith 
822811af0c4SBarry Smith    Notes:
82320f4b53cSBarry Smith    Since `PetscOptionsInsert()` is automatically called by `PetscInitialize()`,
824811af0c4SBarry Smith    the user does not typically need to call this routine. `PetscOptionsInsert()`
825811af0c4SBarry Smith    can be called several times, adding additional entries into the database.
826811af0c4SBarry Smith 
827811af0c4SBarry Smith    See `PetscInitialize()` for options related to option database monitoring.
828081c24baSBoyana Norris 
829db781477SPatrick Sanan .seealso: `PetscOptionsDestroy()`, `PetscOptionsView()`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`,
830db781477SPatrick Sanan           `PetscInitialize()`
831e5c89e4eSSatish Balay @*/
832d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsert(PetscOptions options, int *argc, char ***args, const char file[])
833d71ae5a4SJacob Faibussowitsch {
834d06005cbSLisandro Dalcin   MPI_Comm    comm = PETSC_COMM_WORLD;
835e5c89e4eSSatish Balay   PetscMPIInt rank;
836c5b5d8d5SVaclav Hapla   PetscBool   hasArgs     = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
837c5b5d8d5SVaclav Hapla   PetscBool   skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;
838e5c89e4eSSatish Balay 
839e5c89e4eSSatish Balay   PetscFunctionBegin;
84008401ef6SPierre Jolivet   PetscCheck(!hasArgs || (args && *args), comm, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
8419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
842e5c89e4eSSatish Balay 
843c5b5d8d5SVaclav Hapla   if (!options) {
8449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsCreateDefault());
845c5b5d8d5SVaclav Hapla     options = defaultoptions;
846c5b5d8d5SVaclav Hapla   }
847c5b5d8d5SVaclav Hapla   if (hasArgs) {
848c5b5d8d5SVaclav Hapla     /* process options with absolute precedence */
8499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsProcessPrecedentFlags(options, *argc, *args, &skipPetscrc, &skipPetscrcSet));
850660278c0SBarry Smith     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci", &PetscCIEnabled, NULL));
851c5b5d8d5SVaclav Hapla   }
8524b09e917SBarry Smith   if (file && file[0]) {
8539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, file, PETSC_TRUE));
854c5b5d8d5SVaclav Hapla     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
8559566063dSJacob Faibussowitsch     if (!skipPetscrcSet) PetscCall(PetscOptionsGetBool(options, NULL, "-skip_petscrc", &skipPetscrc, NULL));
856321366bcSBarry Smith   }
857c5b5d8d5SVaclav Hapla   if (!skipPetscrc) {
858be10d61cSLisandro Dalcin     char filename[PETSC_MAX_PATH_LEN];
8599566063dSJacob Faibussowitsch     PetscCall(PetscGetHomeDirectory(filename, sizeof(filename)));
8609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(filename, (int)sizeof(filename), MPI_CHAR, 0, comm));
861c6a7a370SJeremy L Thompson     if (filename[0]) PetscCall(PetscStrlcat(filename, "/.petscrc", sizeof(filename)));
8629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, filename, PETSC_FALSE));
8639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, ".petscrc", PETSC_FALSE));
8649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, "petscrc", PETSC_FALSE));
865e5c89e4eSSatish Balay   }
866e5c89e4eSSatish Balay 
8672d747510SLisandro Dalcin   /* insert environment options */
868e5c89e4eSSatish Balay   {
8692d747510SLisandro Dalcin     char  *eoptions = NULL;
870e5c89e4eSSatish Balay     size_t len      = 0;
871dd400576SPatrick Sanan     if (rank == 0) {
872e5c89e4eSSatish Balay       eoptions = (char *)getenv("PETSC_OPTIONS");
8739566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(eoptions, &len));
874e5c89e4eSSatish Balay     }
8759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
876e5c89e4eSSatish Balay     if (len) {
8779566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
8789566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
87996fc60bcSBarry Smith       if (rank) eoptions[len] = 0;
8809355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertString_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
8819566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscFree(eoptions));
882e5c89e4eSSatish Balay     }
883e5c89e4eSSatish Balay   }
884e5c89e4eSSatish Balay 
885d06005cbSLisandro Dalcin   /* insert YAML environment options */
88656a31166SBarry Smith   {
8879fc438c3SToby Isaac     char  *eoptions = NULL;
8889fc438c3SToby Isaac     size_t len      = 0;
889dd400576SPatrick Sanan     if (rank == 0) {
8909fc438c3SToby Isaac       eoptions = (char *)getenv("PETSC_OPTIONS_YAML");
8919566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(eoptions, &len));
8929fc438c3SToby Isaac     }
8939566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
8949fc438c3SToby Isaac     if (len) {
8959566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
8969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
8979fc438c3SToby Isaac       if (rank) eoptions[len] = 0;
8989355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
8999566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscFree(eoptions));
9009fc438c3SToby Isaac     }
9019fc438c3SToby Isaac   }
9023bcbd388SSean Farley 
903c5b5d8d5SVaclav Hapla   /* insert command line options here because they take precedence over arguments in petscrc/environment */
9049566063dSJacob Faibussowitsch   if (hasArgs) PetscCall(PetscOptionsInsertArgs(options, *argc - 1, *args + 1));
905660278c0SBarry Smith   PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci_portable_error_output", &PetscCIEnabledPortableErrorOutput, NULL));
9063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
907e5c89e4eSSatish Balay }
908e5c89e4eSSatish Balay 
909660278c0SBarry Smith /* These options are not printed with PetscOptionsView() or PetscOptionsMonitor() when PetscCIEnabled is on */
910660278c0SBarry Smith /* TODO: get the list from the test harness, do not have it hardwired here. Maybe from gmakegentest.py */
9119371c9d4SSatish Balay static const char *PetscCIOptions[] = {
9129371c9d4SSatish Balay   "malloc_debug",
913660278c0SBarry Smith   "malloc_dump",
914660278c0SBarry Smith   "malloc_test",
9151f4c7579SBarry Smith   "malloc",
916660278c0SBarry Smith   "nox",
917660278c0SBarry Smith   "nox_warning",
918660278c0SBarry Smith   "display",
919660278c0SBarry Smith   "saws_port_auto_select",
920660278c0SBarry Smith   "saws_port_auto_select_silent",
921660278c0SBarry Smith   "vecscatter_mpi1",
922660278c0SBarry Smith   "check_pointer_intensity",
923660278c0SBarry Smith   "cuda_initialize",
924660278c0SBarry Smith   "error_output_stdout",
925660278c0SBarry Smith   "use_gpu_aware_mpi",
926660278c0SBarry Smith   "checkfunctionlist",
9273872ee93SStefano Zampini   "fp_trap",
928660278c0SBarry Smith   "petsc_ci",
929660278c0SBarry Smith   "petsc_ci_portable_error_output",
930660278c0SBarry Smith };
931660278c0SBarry Smith 
932d71ae5a4SJacob Faibussowitsch static PetscBool PetscCIOption(const char *name)
933d71ae5a4SJacob Faibussowitsch {
934660278c0SBarry Smith   PetscInt  idx;
935660278c0SBarry Smith   PetscBool found;
936660278c0SBarry Smith 
937660278c0SBarry Smith   if (!PetscCIEnabled) return PETSC_FALSE;
9383ba16761SJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscEListFind(PETSC_STATIC_ARRAY_LENGTH(PetscCIOptions), PetscCIOptions, name, &idx, &found));
939660278c0SBarry Smith   return found;
940660278c0SBarry Smith }
941660278c0SBarry Smith 
942e5c89e4eSSatish Balay /*@C
94388c29154SBarry Smith    PetscOptionsView - Prints the options that have been loaded. This is
944e5c89e4eSSatish Balay    useful for debugging purposes.
945e5c89e4eSSatish Balay 
946c3339decSBarry Smith    Logically Collective
947e5c89e4eSSatish Balay 
948d8d19677SJose E. Roman    Input Parameters:
94920f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
950811af0c4SBarry Smith -  viewer - must be an `PETSCVIEWERASCII` viewer
951e5c89e4eSSatish Balay 
952e5c89e4eSSatish Balay    Options Database Key:
953811af0c4SBarry Smith .  -options_view - Activates `PetscOptionsView()` within `PetscFinalize()`
954e5c89e4eSSatish Balay 
95520f4b53cSBarry Smith    Level: advanced
95620f4b53cSBarry Smith 
957811af0c4SBarry Smith    Note:
95821532e8aSBarry Smith    Only the MPI rank 0 of the `MPI_Comm` used to create view prints the option values. Other processes
9591c9f3c13SBarry Smith    may have different values but they are not printed.
9601c9f3c13SBarry Smith 
961db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`
962e5c89e4eSSatish Balay @*/
963d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsView(PetscOptions options, PetscViewer viewer)
964d71ae5a4SJacob Faibussowitsch {
965660278c0SBarry Smith   PetscInt  i, N = 0;
96688c29154SBarry Smith   PetscBool isascii;
967e5c89e4eSSatish Balay 
968e5c89e4eSSatish Balay   PetscFunctionBegin;
9692d747510SLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
970c5929fdfSBarry Smith   options = options ? options : defaultoptions;
97188c29154SBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
9729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
97328b400f6SJacob Faibussowitsch   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Only supports ASCII viewer");
97488c29154SBarry Smith 
975660278c0SBarry Smith   for (i = 0; i < options->N; i++) {
976660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
977660278c0SBarry Smith     N++;
978660278c0SBarry Smith   }
979660278c0SBarry Smith 
980660278c0SBarry Smith   if (!N) {
9819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "#No PETSc Option Table entries\n"));
9823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
98330694fe9SBarry Smith   }
9842d747510SLisandro Dalcin 
9859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "#PETSc Option Table entries:\n"));
986e5c89e4eSSatish Balay   for (i = 0; i < options->N; i++) {
987660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
988e5c89e4eSSatish Balay     if (options->values[i]) {
9899355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s %s", options->names[i], options->values[i]));
990e5c89e4eSSatish Balay     } else {
9919355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s", options->names[i]));
992e5c89e4eSSatish Balay     }
9939355ec05SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, " # (source: %s)\n", PetscOptionSources[options->source[i]]));
994e5c89e4eSSatish Balay   }
9959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "#End of PETSc Option Table entries\n"));
9963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
997e5c89e4eSSatish Balay }
998e5c89e4eSSatish Balay 
999e11779c2SBarry Smith /*
1000e11779c2SBarry Smith    Called by error handlers to print options used in run
1001e11779c2SBarry Smith */
1002d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftError(void)
1003d71ae5a4SJacob Faibussowitsch {
1004f4bc716fSBarry Smith   PetscInt i, nopt = 0;
1005f4bc716fSBarry Smith 
1006f4bc716fSBarry Smith   for (i = 0; i < defaultoptions->N; i++) {
1007f4bc716fSBarry Smith     if (!defaultoptions->used[i]) {
1008f4bc716fSBarry Smith       if (PetscCIOption(defaultoptions->names[i])) continue;
1009f4bc716fSBarry Smith       nopt++;
1010f4bc716fSBarry Smith     }
1011f4bc716fSBarry Smith   }
1012f4bc716fSBarry Smith   if (nopt) {
10137d44c3adSBarry Smith     PetscCall((*PetscErrorPrintf)("WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!\n"));
1014f4bc716fSBarry Smith     for (i = 0; i < defaultoptions->N; i++) {
1015f4bc716fSBarry Smith       if (!defaultoptions->used[i]) {
1016f4bc716fSBarry Smith         if (PetscCIOption(defaultoptions->names[i])) continue;
10173ba16761SJacob Faibussowitsch         if (defaultoptions->values[i]) PetscCall((*PetscErrorPrintf)("  Option left: name:-%s value: %s source: %s\n", defaultoptions->names[i], defaultoptions->values[i], PetscOptionSources[defaultoptions->source[i]]));
10183ba16761SJacob Faibussowitsch         else PetscCall((*PetscErrorPrintf)("  Option left: name:-%s (no value) source: %s\n", defaultoptions->names[i], PetscOptionSources[defaultoptions->source[i]]));
1019f4bc716fSBarry Smith       }
1020f4bc716fSBarry Smith     }
1021f4bc716fSBarry Smith   }
10223ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1023f4bc716fSBarry Smith }
1024f4bc716fSBarry Smith 
1025d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
1026d71ae5a4SJacob Faibussowitsch {
1027660278c0SBarry Smith   PetscInt     i, N = 0;
10284416b707SBarry Smith   PetscOptions options = defaultoptions;
1029e11779c2SBarry Smith 
1030660278c0SBarry Smith   for (i = 0; i < options->N; i++) {
1031660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
1032660278c0SBarry Smith     N++;
1033660278c0SBarry Smith   }
1034660278c0SBarry Smith 
1035660278c0SBarry Smith   if (N) {
10363ba16761SJacob Faibussowitsch     PetscCall((*PetscErrorPrintf)("PETSc Option Table entries:\n"));
1037e11779c2SBarry Smith   } else {
10383ba16761SJacob Faibussowitsch     PetscCall((*PetscErrorPrintf)("No PETSc Option Table entries\n"));
1039e11779c2SBarry Smith   }
1040e11779c2SBarry Smith   for (i = 0; i < options->N; i++) {
1041660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
1042e11779c2SBarry Smith     if (options->values[i]) {
10433ba16761SJacob Faibussowitsch       PetscCall((*PetscErrorPrintf)("-%s %s (source: %s)\n", options->names[i], options->values[i], PetscOptionSources[options->source[i]]));
1044e11779c2SBarry Smith     } else {
10453ba16761SJacob Faibussowitsch       PetscCall((*PetscErrorPrintf)("-%s (source: %s)\n", options->names[i], PetscOptionSources[options->source[i]]));
1046e11779c2SBarry Smith     }
1047e11779c2SBarry Smith   }
10483ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1049e11779c2SBarry Smith }
1050e11779c2SBarry Smith 
1051e5c89e4eSSatish Balay /*@C
105274e0666dSJed Brown    PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.
105374e0666dSJed Brown 
10541c9f3c13SBarry Smith    Logically Collective
105574e0666dSJed Brown 
1056d8d19677SJose E. Roman    Input Parameters:
105720f4b53cSBarry Smith +  options - options database, or `NULL` for the default global database
1058c5929fdfSBarry Smith -  prefix - The string to append to the existing prefix
10599db968c8SJed Brown 
10609db968c8SJed Brown    Options Database Keys:
10619db968c8SJed Brown +   -prefix_push <some_prefix_> - push the given prefix
10629db968c8SJed Brown -   -prefix_pop - pop the last prefix
10639db968c8SJed Brown 
106420f4b53cSBarry Smith    Level: advanced
106520f4b53cSBarry Smith 
10669db968c8SJed Brown    Notes:
106721532e8aSBarry Smith    It is common to use this in conjunction with `-options_file` as in
10689db968c8SJed Brown 
10699db968c8SJed Brown $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
10709db968c8SJed Brown 
107121532e8aSBarry Smith    where the files no longer require all options to be prefixed with `-system2_`.
107274e0666dSJed Brown 
107321532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes that call this routine will
10741c9f3c13SBarry Smith    have the affect of these options. If some processes that create objects call this routine and others do
10751c9f3c13SBarry Smith    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
10761c9f3c13SBarry Smith    on different ranks.
10771c9f3c13SBarry Smith 
1078db781477SPatrick Sanan .seealso: `PetscOptionsPrefixPop()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
107974e0666dSJed Brown @*/
1080d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPrefixPush(PetscOptions options, const char prefix[])
1081d71ae5a4SJacob Faibussowitsch {
108274e0666dSJed Brown   size_t    n;
108374e0666dSJed Brown   PetscInt  start;
10849355ec05SMatthew G. Knepley   char      key[PETSC_MAX_OPTION_NAME + 1];
10852d747510SLisandro Dalcin   PetscBool valid;
108674e0666dSJed Brown 
108774e0666dSJed Brown   PetscFunctionBegin;
1088064a246eSJacob Faibussowitsch   PetscValidCharPointer(prefix, 2);
1089c5929fdfSBarry Smith   options = options ? options : defaultoptions;
109008401ef6SPierre Jolivet   PetscCheck(options->prefixind < MAXPREFIXES, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum depth of prefix stack %d exceeded, recompile \n src/sys/objects/options.c with larger value for MAXPREFIXES", MAXPREFIXES);
10912d747510SLisandro Dalcin   key[0] = '-'; /* keys must start with '-' */
10929566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(key + 1, prefix, sizeof(key) - 1));
10939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(key, &valid));
10948bf569ecSLisandro Dalcin   if (!valid && options->prefixind > 0 && isdigit((int)prefix[0])) valid = PETSC_TRUE; /* If the prefix stack is not empty, make numbers a valid prefix */
109528b400f6SJacob Faibussowitsch   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_USER, "Given prefix \"%s\" not valid (the first character must be a letter%s, do not include leading '-')", prefix, options->prefixind ? " or digit" : "");
109674e0666dSJed Brown   start = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
10979566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(prefix, &n));
109808401ef6SPierre Jolivet   PetscCheck(n + 1 <= sizeof(options->prefix) - start, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum prefix length %zu exceeded", sizeof(options->prefix));
10999566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(options->prefix + start, prefix, n + 1));
110074e0666dSJed Brown   options->prefixstack[options->prefixind++] = start + n;
11013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110274e0666dSJed Brown }
110374e0666dSJed Brown 
1104c5929fdfSBarry Smith /*@C
1105811af0c4SBarry Smith    PetscOptionsPrefixPop - Remove the latest options prefix, see `PetscOptionsPrefixPush()` for details
110674e0666dSJed Brown 
1107811af0c4SBarry Smith    Logically Collective on the `MPI_Comm` used when called `PetscOptionsPrefixPush()`
110874e0666dSJed Brown 
1109811af0c4SBarry Smith   Input Parameter:
111020f4b53cSBarry Smith .  options - options database, or `NULL` for the default global database
1111c5929fdfSBarry Smith 
111274e0666dSJed Brown    Level: advanced
111374e0666dSJed Brown 
1114db781477SPatrick Sanan .seealso: `PetscOptionsPrefixPush()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
111574e0666dSJed Brown @*/
1116d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1117d71ae5a4SJacob Faibussowitsch {
111874e0666dSJed Brown   PetscInt offset;
111974e0666dSJed Brown 
112074e0666dSJed Brown   PetscFunctionBegin;
1121c5929fdfSBarry Smith   options = options ? options : defaultoptions;
112208401ef6SPierre Jolivet   PetscCheck(options->prefixind >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "More prefixes popped than pushed");
112374e0666dSJed Brown   options->prefixind--;
112474e0666dSJed Brown   offset                  = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
112574e0666dSJed Brown   options->prefix[offset] = 0;
11263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112774e0666dSJed Brown }
112874e0666dSJed Brown 
1129a542b6e8SBarry Smith /*@C
1130a542b6e8SBarry Smith     PetscOptionsClear - Removes all options form the database leaving it empty.
1131a542b6e8SBarry Smith 
11321c9f3c13SBarry Smith     Logically Collective
11331c9f3c13SBarry Smith 
1134811af0c4SBarry Smith   Input Parameter:
113520f4b53cSBarry Smith .  options - options database, use `NULL` for the default global database
1136c5929fdfSBarry Smith 
113720f4b53cSBarry Smith    Level: developer
113820f4b53cSBarry Smith 
113920f4b53cSBarry Smith    Note:
114021532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes that call this routine will
11411c9f3c13SBarry Smith    have the affect of these options. If some processes that create objects call this routine and others do
11421c9f3c13SBarry Smith    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
11431c9f3c13SBarry Smith    on different ranks.
11441c9f3c13SBarry Smith 
1145db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`
1146a542b6e8SBarry Smith @*/
1147d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsClear(PetscOptions options)
1148d71ae5a4SJacob Faibussowitsch {
1149a542b6e8SBarry Smith   PetscInt i;
1150a542b6e8SBarry Smith 
115139a651e2SJacob Faibussowitsch   PetscFunctionBegin;
1152c5929fdfSBarry Smith   options = options ? options : defaultoptions;
11533ba16761SJacob Faibussowitsch   if (!options) PetscFunctionReturn(PETSC_SUCCESS);
11542d747510SLisandro Dalcin 
1155a542b6e8SBarry Smith   for (i = 0; i < options->N; i++) {
1156a542b6e8SBarry Smith     if (options->names[i]) free(options->names[i]);
1157a542b6e8SBarry Smith     if (options->values[i]) free(options->values[i]);
1158a542b6e8SBarry Smith   }
11592d747510SLisandro Dalcin   options->N = 0;
11609355ec05SMatthew G. Knepley   free(options->names);
11619355ec05SMatthew G. Knepley   free(options->values);
11629355ec05SMatthew G. Knepley   free(options->used);
11639355ec05SMatthew G. Knepley   free(options->source);
11649355ec05SMatthew G. Knepley   options->names  = NULL;
11659355ec05SMatthew G. Knepley   options->values = NULL;
11669355ec05SMatthew G. Knepley   options->used   = NULL;
11679355ec05SMatthew G. Knepley   options->source = NULL;
11689355ec05SMatthew G. Knepley   options->Nalloc = 0;
11692d747510SLisandro Dalcin 
11709355ec05SMatthew G. Knepley   for (i = 0; i < options->Na; i++) {
1171a542b6e8SBarry Smith     free(options->aliases1[i]);
1172a542b6e8SBarry Smith     free(options->aliases2[i]);
1173a542b6e8SBarry Smith   }
11749355ec05SMatthew G. Knepley   options->Na = 0;
11759355ec05SMatthew G. Knepley   free(options->aliases1);
11769355ec05SMatthew G. Knepley   free(options->aliases2);
11779355ec05SMatthew G. Knepley   options->aliases1 = options->aliases2 = NULL;
11789355ec05SMatthew G. Knepley   options->Naalloc                      = 0;
1179a542b6e8SBarry Smith 
11802d747510SLisandro Dalcin   /* destroy hash table */
11812d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
11822d747510SLisandro Dalcin   options->ht = NULL;
11830eb63584SBarry Smith 
11842d747510SLisandro Dalcin   options->prefixind  = 0;
11852d747510SLisandro Dalcin   options->prefix[0]  = 0;
11862d747510SLisandro Dalcin   options->help       = PETSC_FALSE;
11879355ec05SMatthew G. Knepley   options->help_intro = PETSC_FALSE;
11883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11894416b707SBarry Smith }
11904416b707SBarry Smith 
11912d747510SLisandro Dalcin /*@C
11922d747510SLisandro Dalcin    PetscOptionsSetAlias - Makes a key and alias for another key
11932d747510SLisandro Dalcin 
11941c9f3c13SBarry Smith    Logically Collective
11952d747510SLisandro Dalcin 
11962d747510SLisandro Dalcin    Input Parameters:
119720f4b53cSBarry Smith +  options - options database, or `NULL` for default global database
11982d747510SLisandro Dalcin .  newname - the alias
11992d747510SLisandro Dalcin -  oldname - the name that alias will refer to
12002d747510SLisandro Dalcin 
12012d747510SLisandro Dalcin    Level: advanced
12022d747510SLisandro Dalcin 
120320f4b53cSBarry Smith    Note:
120421532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes that call this routine will
12051c9f3c13SBarry Smith    have the affect of these options. If some processes that create objects call this routine and others do
12061c9f3c13SBarry Smith    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
12071c9f3c13SBarry Smith    on different ranks.
12081c9f3c13SBarry Smith 
1209c2e3fba1SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1210c2e3fba1SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1211db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1212c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1213db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1214db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
12152d747510SLisandro Dalcin @*/
1216d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSetAlias(PetscOptions options, const char newname[], const char oldname[])
1217d71ae5a4SJacob Faibussowitsch {
12182d747510SLisandro Dalcin   size_t    len;
12199210b8eaSVaclav Hapla   PetscBool valid;
12202d747510SLisandro Dalcin 
12212d747510SLisandro Dalcin   PetscFunctionBegin;
12222d747510SLisandro Dalcin   PetscValidCharPointer(newname, 2);
12232d747510SLisandro Dalcin   PetscValidCharPointer(oldname, 3);
12242d747510SLisandro Dalcin   options = options ? options : defaultoptions;
12259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(newname, &valid));
122628b400f6SJacob Faibussowitsch   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliased option %s", newname);
12279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(oldname, &valid));
122828b400f6SJacob Faibussowitsch   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliasee option %s", oldname);
12292d747510SLisandro Dalcin 
12309355ec05SMatthew G. Knepley   if (options->Na == options->Naalloc) {
12319355ec05SMatthew G. Knepley     char **tmpA1, **tmpA2;
12322d747510SLisandro Dalcin 
12339355ec05SMatthew G. Knepley     options->Naalloc = PetscMax(4, options->Naalloc * 2);
12349355ec05SMatthew G. Knepley     tmpA1            = (char **)malloc(options->Naalloc * sizeof(char *));
12359355ec05SMatthew G. Knepley     tmpA2            = (char **)malloc(options->Naalloc * sizeof(char *));
12369355ec05SMatthew G. Knepley     for (int i = 0; i < options->Na; ++i) {
12379355ec05SMatthew G. Knepley       tmpA1[i] = options->aliases1[i];
12389355ec05SMatthew G. Knepley       tmpA2[i] = options->aliases2[i];
12399355ec05SMatthew G. Knepley     }
12409355ec05SMatthew G. Knepley     free(options->aliases1);
12419355ec05SMatthew G. Knepley     free(options->aliases2);
12429355ec05SMatthew G. Knepley     options->aliases1 = tmpA1;
12439355ec05SMatthew G. Knepley     options->aliases2 = tmpA2;
12449355ec05SMatthew G. Knepley   }
12459371c9d4SSatish Balay   newname++;
12469371c9d4SSatish Balay   oldname++;
12479566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(newname, &len));
12489355ec05SMatthew G. Knepley   options->aliases1[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1249c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(options->aliases1[options->Na], newname, len + 1));
12509566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(oldname, &len));
12519355ec05SMatthew G. Knepley   options->aliases2[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1252c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(options->aliases2[options->Na], oldname, len + 1));
12539355ec05SMatthew G. Knepley   ++options->Na;
12543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12552d747510SLisandro Dalcin }
12564416b707SBarry Smith 
1257e5c89e4eSSatish Balay /*@C
1258e5c89e4eSSatish Balay    PetscOptionsSetValue - Sets an option name-value pair in the options
1259e5c89e4eSSatish Balay    database, overriding whatever is already present.
1260e5c89e4eSSatish Balay 
12611c9f3c13SBarry Smith    Logically Collective
1262e5c89e4eSSatish Balay 
1263e5c89e4eSSatish Balay    Input Parameters:
126420f4b53cSBarry Smith +  options - options database, use `NULL` for the default global database
1265c5929fdfSBarry Smith .  name - name of option, this SHOULD have the - prepended
126620f4b53cSBarry Smith -  value - the option value (not used for all options, so can be `NULL`)
1267e5c89e4eSSatish Balay 
1268e5c89e4eSSatish Balay    Level: intermediate
1269e5c89e4eSSatish Balay 
1270e5c89e4eSSatish Balay    Note:
1271811af0c4SBarry Smith    This function can be called BEFORE `PetscInitialize()`
1272d49172ceSBarry Smith 
127321532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes that call this routine will
12741c9f3c13SBarry Smith    have the affect of these options. If some processes that create objects call this routine and others do
12751c9f3c13SBarry Smith    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
12761c9f3c13SBarry Smith    on different ranks.
12771c9f3c13SBarry Smith 
127820f4b53cSBarry Smith    Developers Note:
127920f4b53cSBarry Smith    Uses malloc() directly because PETSc may not be initialized yet.
1280b0250c70SBarry Smith 
1281db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`, `PetscOptionsClearValue()`
1282e5c89e4eSSatish Balay @*/
1283d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSetValue(PetscOptions options, const char name[], const char value[])
1284d71ae5a4SJacob Faibussowitsch {
128539a651e2SJacob Faibussowitsch   PetscFunctionBegin;
12869355ec05SMatthew G. Knepley   PetscCall(PetscOptionsSetValue_Private(options, name, value, NULL, PETSC_OPT_CODE));
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288c5b5d8d5SVaclav Hapla }
1289c5b5d8d5SVaclav Hapla 
12909355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options, const char name[], const char value[], int *pos, PetscOptionSource source)
1291d71ae5a4SJacob Faibussowitsch {
1292e5c89e4eSSatish Balay   size_t    len;
12939355ec05SMatthew G. Knepley   int       n, i;
1294e5c89e4eSSatish Balay   char    **names;
12959355ec05SMatthew G. Knepley   char      fullname[PETSC_MAX_OPTION_NAME] = "";
1296c5b5d8d5SVaclav Hapla   PetscBool flg;
1297e5c89e4eSSatish Balay 
129839a651e2SJacob Faibussowitsch   PetscFunctionBegin;
12997272c0d2SVaclav Hapla   if (!options) {
13009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsCreateDefault());
13017272c0d2SVaclav Hapla     options = defaultoptions;
1302c5929fdfSBarry Smith   }
130339a651e2SJacob Faibussowitsch   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "name %s must start with '-'", name);
13042d747510SLisandro Dalcin 
13059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsSkipPrecedent(options, name, &flg));
13063ba16761SJacob Faibussowitsch   if (flg) PetscFunctionReturn(PETSC_SUCCESS);
1307e5c89e4eSSatish Balay 
13082d747510SLisandro Dalcin   name++; /* skip starting dash */
13092d747510SLisandro Dalcin 
131074e0666dSJed Brown   if (options->prefixind > 0) {
1311d49172ceSBarry Smith     strncpy(fullname, options->prefix, sizeof(fullname));
13122d747510SLisandro Dalcin     fullname[sizeof(fullname) - 1] = 0;
131389ae1891SBarry Smith     strncat(fullname, name, sizeof(fullname) - strlen(fullname) - 1);
13142d747510SLisandro Dalcin     fullname[sizeof(fullname) - 1] = 0;
131574e0666dSJed Brown     name                           = fullname;
131674e0666dSJed Brown   }
131774e0666dSJed Brown 
131874e0666dSJed Brown   /* check against aliases */
13199355ec05SMatthew G. Knepley   for (i = 0; i < options->Na; i++) {
13202d747510SLisandro Dalcin     int result = PetscOptNameCmp(options->aliases1[i], name);
13219371c9d4SSatish Balay     if (!result) {
13229371c9d4SSatish Balay       name = options->aliases2[i];
13239371c9d4SSatish Balay       break;
13249371c9d4SSatish Balay     }
1325e5c89e4eSSatish Balay   }
1326e5c89e4eSSatish Balay 
13272d747510SLisandro Dalcin   /* slow search */
13289355ec05SMatthew G. Knepley   n     = options->N;
1329e5c89e4eSSatish Balay   names = options->names;
13309355ec05SMatthew G. Knepley   for (i = 0; i < options->N; i++) {
13312d747510SLisandro Dalcin     int result = PetscOptNameCmp(names[i], name);
13322d747510SLisandro Dalcin     if (!result) {
13339371c9d4SSatish Balay       n = i;
13349371c9d4SSatish Balay       goto setvalue;
13352d747510SLisandro Dalcin     } else if (result > 0) {
13369371c9d4SSatish Balay       n = i;
13379371c9d4SSatish Balay       break;
1338e5c89e4eSSatish Balay     }
1339e5c89e4eSSatish Balay   }
13409355ec05SMatthew G. Knepley   if (options->N == options->Nalloc) {
13419355ec05SMatthew G. Knepley     char             **names, **values;
13429355ec05SMatthew G. Knepley     PetscBool         *used;
13439355ec05SMatthew G. Knepley     PetscOptionSource *source;
13449355ec05SMatthew G. Knepley 
13459355ec05SMatthew G. Knepley     options->Nalloc = PetscMax(10, options->Nalloc * 2);
13469355ec05SMatthew G. Knepley     names           = (char **)malloc(options->Nalloc * sizeof(char *));
13479355ec05SMatthew G. Knepley     values          = (char **)malloc(options->Nalloc * sizeof(char *));
13489355ec05SMatthew G. Knepley     used            = (PetscBool *)malloc(options->Nalloc * sizeof(PetscBool));
13499355ec05SMatthew G. Knepley     source          = (PetscOptionSource *)malloc(options->Nalloc * sizeof(PetscOptionSource));
13509355ec05SMatthew G. Knepley     for (int i = 0; i < options->N; ++i) {
13519355ec05SMatthew G. Knepley       names[i]  = options->names[i];
13529355ec05SMatthew G. Knepley       values[i] = options->values[i];
13539355ec05SMatthew G. Knepley       used[i]   = options->used[i];
13549355ec05SMatthew G. Knepley       source[i] = options->source[i];
13559355ec05SMatthew G. Knepley     }
13569355ec05SMatthew G. Knepley     free(options->names);
13579355ec05SMatthew G. Knepley     free(options->values);
13589355ec05SMatthew G. Knepley     free(options->used);
13599355ec05SMatthew G. Knepley     free(options->source);
13609355ec05SMatthew G. Knepley     options->names  = names;
13619355ec05SMatthew G. Knepley     options->values = values;
13629355ec05SMatthew G. Knepley     options->used   = used;
13639355ec05SMatthew G. Knepley     options->source = source;
13649355ec05SMatthew G. Knepley   }
136539a651e2SJacob Faibussowitsch 
13662d747510SLisandro Dalcin   /* shift remaining values up 1 */
13679355ec05SMatthew G. Knepley   for (i = options->N; i > n; i--) {
13685e8c5e88SLisandro Dalcin     options->names[i]  = options->names[i - 1];
1369e5c89e4eSSatish Balay     options->values[i] = options->values[i - 1];
1370e5c89e4eSSatish Balay     options->used[i]   = options->used[i - 1];
13719355ec05SMatthew G. Knepley     options->source[i] = options->source[i - 1];
1372e5c89e4eSSatish Balay   }
13732d747510SLisandro Dalcin   options->names[n]  = NULL;
13742d747510SLisandro Dalcin   options->values[n] = NULL;
13752d747510SLisandro Dalcin   options->used[n]   = PETSC_FALSE;
13769355ec05SMatthew G. Knepley   options->source[n] = PETSC_OPT_CODE;
13772d747510SLisandro Dalcin   options->N++;
13782d747510SLisandro Dalcin 
13792d747510SLisandro Dalcin   /* destroy hash table */
13802d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
13812d747510SLisandro Dalcin   options->ht = NULL;
13822d747510SLisandro Dalcin 
13832d747510SLisandro Dalcin   /* set new name */
138470d8d27cSBarry Smith   len               = strlen(name);
13855e8c5e88SLisandro Dalcin   options->names[n] = (char *)malloc((len + 1) * sizeof(char));
138639a651e2SJacob Faibussowitsch   PetscCheck(options->names[n], PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate option name");
1387d49172ceSBarry Smith   strcpy(options->names[n], name);
13882d747510SLisandro Dalcin 
13892d747510SLisandro Dalcin setvalue:
13902d747510SLisandro Dalcin   /* set new value */
13912d747510SLisandro Dalcin   if (options->values[n]) free(options->values[n]);
1392d49172ceSBarry Smith   len = value ? strlen(value) : 0;
13935e8c5e88SLisandro Dalcin   if (len) {
1394e5c89e4eSSatish Balay     options->values[n] = (char *)malloc((len + 1) * sizeof(char));
1395d49172ceSBarry Smith     if (!options->values[n]) return PETSC_ERR_MEM;
1396d49172ceSBarry Smith     strcpy(options->values[n], value);
13972d747510SLisandro Dalcin   } else {
13982d747510SLisandro Dalcin     options->values[n] = NULL;
13992d747510SLisandro Dalcin   }
14009355ec05SMatthew G. Knepley   options->source[n] = source;
14012d747510SLisandro Dalcin 
140291ad3481SVaclav Hapla   /* handle -help so that it can be set from anywhere */
140391ad3481SVaclav Hapla   if (!PetscOptNameCmp(name, "help")) {
140491ad3481SVaclav Hapla     options->help       = PETSC_TRUE;
1405d06005cbSLisandro Dalcin     options->help_intro = (value && !PetscOptNameCmp(value, "intro")) ? PETSC_TRUE : PETSC_FALSE;
140691ad3481SVaclav Hapla     options->used[n]    = PETSC_TRUE;
140791ad3481SVaclav Hapla   }
140891ad3481SVaclav Hapla 
14099355ec05SMatthew G. Knepley   PetscCall(PetscOptionsMonitor(options, name, value, source));
1410c5b5d8d5SVaclav Hapla   if (pos) *pos = n;
14113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1412e5c89e4eSSatish Balay }
1413e5c89e4eSSatish Balay 
1414e5c89e4eSSatish Balay /*@C
1415e5c89e4eSSatish Balay    PetscOptionsClearValue - Clears an option name-value pair in the options
1416e5c89e4eSSatish Balay    database, overriding whatever is already present.
1417e5c89e4eSSatish Balay 
14181c9f3c13SBarry Smith    Logically Collective
1419e5c89e4eSSatish Balay 
1420d8d19677SJose E. Roman    Input Parameters:
142120f4b53cSBarry Smith +  options - options database, use `NULL` for the default global database
1422a2b725a8SWilliam Gropp -  name - name of option, this SHOULD have the - prepended
1423e5c89e4eSSatish Balay 
1424e5c89e4eSSatish Balay    Level: intermediate
1425e5c89e4eSSatish Balay 
1426811af0c4SBarry Smith    Note:
142721532e8aSBarry Smith    The collectivity of this routine is complex; only the MPI processes that call this routine will
14281c9f3c13SBarry Smith    have the affect of these options. If some processes that create objects call this routine and others do
14291c9f3c13SBarry Smith    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
14301c9f3c13SBarry Smith    on different ranks.
14311c9f3c13SBarry Smith 
1432db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`
1433e5c89e4eSSatish Balay @*/
1434d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsClearValue(PetscOptions options, const char name[])
1435d71ae5a4SJacob Faibussowitsch {
14362d747510SLisandro Dalcin   int    N, n, i;
14372d747510SLisandro Dalcin   char **names;
1438e5c89e4eSSatish Balay 
1439e5c89e4eSSatish Balay   PetscFunctionBegin;
1440c5929fdfSBarry Smith   options = options ? options : defaultoptions;
1441cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1442c9dcd962SLisandro Dalcin   if (!PetscOptNameCmp(name, "-help")) options->help = options->help_intro = PETSC_FALSE;
14432d747510SLisandro Dalcin 
14442d747510SLisandro Dalcin   name++; /* skip starting dash */
14452d747510SLisandro Dalcin 
14462d747510SLisandro Dalcin   /* slow search */
14472d747510SLisandro Dalcin   N = n = options->N;
1448e5c89e4eSSatish Balay   names = options->names;
1449e5c89e4eSSatish Balay   for (i = 0; i < N; i++) {
14502d747510SLisandro Dalcin     int result = PetscOptNameCmp(names[i], name);
14512d747510SLisandro Dalcin     if (!result) {
14529371c9d4SSatish Balay       n = i;
14539371c9d4SSatish Balay       break;
14542d747510SLisandro Dalcin     } else if (result > 0) {
14559371c9d4SSatish Balay       n = N;
14569371c9d4SSatish Balay       break;
1457e5c89e4eSSatish Balay     }
14582d747510SLisandro Dalcin   }
14593ba16761SJacob Faibussowitsch   if (n == N) PetscFunctionReturn(PETSC_SUCCESS); /* it was not present */
1460e5c89e4eSSatish Balay 
14612d747510SLisandro Dalcin   /* remove name and value */
14622d747510SLisandro Dalcin   if (options->names[n]) free(options->names[n]);
14632d747510SLisandro Dalcin   if (options->values[n]) free(options->values[n]);
1464e5c89e4eSSatish Balay   /* shift remaining values down 1 */
1465e5c89e4eSSatish Balay   for (i = n; i < N - 1; i++) {
14665e8c5e88SLisandro Dalcin     options->names[i]  = options->names[i + 1];
1467e5c89e4eSSatish Balay     options->values[i] = options->values[i + 1];
1468e5c89e4eSSatish Balay     options->used[i]   = options->used[i + 1];
14699355ec05SMatthew G. Knepley     options->source[i] = options->source[i + 1];
1470e5c89e4eSSatish Balay   }
1471e5c89e4eSSatish Balay   options->N--;
14722d747510SLisandro Dalcin 
14732d747510SLisandro Dalcin   /* destroy hash table */
14742d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
14752d747510SLisandro Dalcin   options->ht = NULL;
14762d747510SLisandro Dalcin 
14779355ec05SMatthew G. Knepley   PetscCall(PetscOptionsMonitor(options, name, NULL, PETSC_OPT_CODE));
14783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1479e5c89e4eSSatish Balay }
1480e5c89e4eSSatish Balay 
1481e5c89e4eSSatish Balay /*@C
14822d747510SLisandro Dalcin    PetscOptionsFindPair - Gets an option name-value pair from the options database.
1483e5c89e4eSSatish Balay 
14842d747510SLisandro Dalcin    Not Collective
1485e5c89e4eSSatish Balay 
1486e5c89e4eSSatish Balay    Input Parameters:
148720f4b53cSBarry Smith +  options - options database, use `NULL` for the default global database
148820f4b53cSBarry Smith .  pre - the string to prepend to the name or `NULL`, this SHOULD NOT have the "-" prepended
14892d747510SLisandro Dalcin -  name - name of option, this SHOULD have the "-" prepended
1490e5c89e4eSSatish Balay 
14912d747510SLisandro Dalcin    Output Parameters:
14922d747510SLisandro Dalcin +  value - the option value (optional, not used for all options)
14932d747510SLisandro Dalcin -  set - whether the option is set (optional)
1494e5c89e4eSSatish Balay 
149520f4b53cSBarry Smith    Level: developer
149620f4b53cSBarry Smith 
1497811af0c4SBarry Smith    Note:
14989666a313SBarry Smith    Each process may find different values or no value depending on how options were inserted into the database
14991c9f3c13SBarry Smith 
1500db781477SPatrick Sanan .seealso: `PetscOptionsSetValue()`, `PetscOptionsClearValue()`
1501e5c89e4eSSatish Balay @*/
1502d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsFindPair(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1503d71ae5a4SJacob Faibussowitsch {
15049355ec05SMatthew G. Knepley   char      buf[PETSC_MAX_OPTION_NAME];
1505daabea38SBarry Smith   PetscBool usehashtable = PETSC_TRUE;
15062d747510SLisandro Dalcin   PetscBool matchnumbers = PETSC_TRUE;
1507e5c89e4eSSatish Balay 
1508e5c89e4eSSatish Balay   PetscFunctionBegin;
1509c5929fdfSBarry Smith   options = options ? options : defaultoptions;
151008401ef6SPierre Jolivet   PetscCheck(!pre || !PetscUnlikely(pre[0] == '-'), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1511cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1512e5c89e4eSSatish Balay 
15132d747510SLisandro Dalcin   name++; /* skip starting dash */
1514e5c89e4eSSatish Balay 
15157cd08cecSJed Brown   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
15162d747510SLisandro Dalcin   if (pre && pre[0]) {
15172d747510SLisandro Dalcin     char *ptr = buf;
15189371c9d4SSatish Balay     if (name[0] == '-') {
15199371c9d4SSatish Balay       *ptr++ = '-';
15209371c9d4SSatish Balay       name++;
15219371c9d4SSatish Balay     }
15229566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(ptr, pre, buf + sizeof(buf) - ptr));
15239566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
15242d747510SLisandro Dalcin     name = buf;
15257cd08cecSJed Brown   }
15262d747510SLisandro Dalcin 
152776bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
15282f828895SJed Brown     PetscBool valid;
15299355ec05SMatthew G. Knepley     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
15309566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
15319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(key, &valid));
153228b400f6SJacob Faibussowitsch     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
15332f828895SJed Brown   }
1534e5c89e4eSSatish Balay 
15352d747510SLisandro Dalcin   if (!options->ht && usehashtable) {
15362d747510SLisandro Dalcin     int          i, ret;
15372d747510SLisandro Dalcin     khiter_t     it;
15382d747510SLisandro Dalcin     khash_t(HO) *ht;
15392d747510SLisandro Dalcin     ht = kh_init(HO);
154028b400f6SJacob Faibussowitsch     PetscCheck(ht, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15412d747510SLisandro Dalcin     ret = kh_resize(HO, ht, options->N * 2); /* twice the required size to reduce risk of collisions */
154228b400f6SJacob Faibussowitsch     PetscCheck(!ret, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15432d747510SLisandro Dalcin     for (i = 0; i < options->N; i++) {
15442d747510SLisandro Dalcin       it = kh_put(HO, ht, options->names[i], &ret);
154508401ef6SPierre Jolivet       PetscCheck(ret == 1, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15462d747510SLisandro Dalcin       kh_val(ht, it) = i;
15472d747510SLisandro Dalcin     }
15482d747510SLisandro Dalcin     options->ht = ht;
15492d747510SLisandro Dalcin   }
15502d747510SLisandro Dalcin 
15519371c9d4SSatish Balay   if (usehashtable) { /* fast search */
15522d747510SLisandro Dalcin     khash_t(HO) *ht = options->ht;
15532d747510SLisandro Dalcin     khiter_t     it = kh_get(HO, ht, name);
15542d747510SLisandro Dalcin     if (it != kh_end(ht)) {
15552d747510SLisandro Dalcin       int i            = kh_val(ht, it);
1556e5c89e4eSSatish Balay       options->used[i] = PETSC_TRUE;
15572d747510SLisandro Dalcin       if (value) *value = options->values[i];
15582d747510SLisandro Dalcin       if (set) *set = PETSC_TRUE;
15593ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
15602d747510SLisandro Dalcin     }
15619371c9d4SSatish Balay   } else { /* slow search */
15622d747510SLisandro Dalcin     int i, N = options->N;
15632d747510SLisandro Dalcin     for (i = 0; i < N; i++) {
1564daabea38SBarry Smith       int result = PetscOptNameCmp(options->names[i], name);
15652d747510SLisandro Dalcin       if (!result) {
15662d747510SLisandro Dalcin         options->used[i] = PETSC_TRUE;
15672d747510SLisandro Dalcin         if (value) *value = options->values[i];
15682d747510SLisandro Dalcin         if (set) *set = PETSC_TRUE;
15693ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
15702d747510SLisandro Dalcin       } else if (result > 0) {
1571e5c89e4eSSatish Balay         break;
1572e5c89e4eSSatish Balay       }
1573e5c89e4eSSatish Balay     }
15742d747510SLisandro Dalcin   }
15752d747510SLisandro Dalcin 
15762d747510SLisandro Dalcin   /*
15772d747510SLisandro Dalcin    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
15782d747510SLisandro Dalcin    Maybe this special lookup mode should be enabled on request with a push/pop API.
15792d747510SLisandro Dalcin    The feature of matching _%d_ used sparingly in the codebase.
15802d747510SLisandro Dalcin    */
15812d747510SLisandro Dalcin   if (matchnumbers) {
15822d747510SLisandro Dalcin     int i, j, cnt = 0, locs[16], loce[16];
1583e5c89e4eSSatish Balay     /* determine the location and number of all _%d_ in the key */
15842d747510SLisandro Dalcin     for (i = 0; name[i]; i++) {
15852d747510SLisandro Dalcin       if (name[i] == '_') {
15862d747510SLisandro Dalcin         for (j = i + 1; name[j]; j++) {
15872d747510SLisandro Dalcin           if (name[j] >= '0' && name[j] <= '9') continue;
15882d747510SLisandro Dalcin           if (name[j] == '_' && j > i + 1) { /* found a number */
1589e5c89e4eSSatish Balay             locs[cnt]   = i + 1;
1590e5c89e4eSSatish Balay             loce[cnt++] = j + 1;
1591e5c89e4eSSatish Balay           }
15922d747510SLisandro Dalcin           i = j - 1;
1593e5c89e4eSSatish Balay           break;
1594e5c89e4eSSatish Balay         }
1595e5c89e4eSSatish Balay       }
1596e5c89e4eSSatish Balay     }
1597e5c89e4eSSatish Balay     for (i = 0; i < cnt; i++) {
15982d747510SLisandro Dalcin       PetscBool found;
15999355ec05SMatthew G. Knepley       char      opt[PETSC_MAX_OPTION_NAME + 1] = "-", tmp[PETSC_MAX_OPTION_NAME];
16009566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(tmp, name, PetscMin((size_t)(locs[i] + 1), sizeof(tmp))));
16019566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, tmp, sizeof(opt)));
16029566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, name + loce[i], sizeof(opt)));
16039566063dSJacob Faibussowitsch       PetscCall(PetscOptionsFindPair(options, NULL, opt, value, &found));
16049371c9d4SSatish Balay       if (found) {
16059371c9d4SSatish Balay         if (set) *set = PETSC_TRUE;
16063ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
16079371c9d4SSatish Balay       }
1608e5c89e4eSSatish Balay     }
1609e5c89e4eSSatish Balay   }
16102d747510SLisandro Dalcin 
16112d747510SLisandro Dalcin   if (set) *set = PETSC_FALSE;
16123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1613e5c89e4eSSatish Balay }
1614e5c89e4eSSatish Balay 
1615d6ced9c0SMatthew G. Knepley /* Check whether any option begins with pre+name */
1616d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1617d71ae5a4SJacob Faibussowitsch {
16189355ec05SMatthew G. Knepley   char buf[PETSC_MAX_OPTION_NAME];
1619d6ced9c0SMatthew G. Knepley   int  numCnt = 0, locs[16], loce[16];
1620514bf10dSMatthew G Knepley 
1621514bf10dSMatthew G Knepley   PetscFunctionBegin;
1622c5929fdfSBarry Smith   options = options ? options : defaultoptions;
1623cc73adaaSBarry Smith   PetscCheck(!pre || pre[0] != '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1624cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1625514bf10dSMatthew G Knepley 
16262d747510SLisandro Dalcin   name++; /* skip starting dash */
1627514bf10dSMatthew G Knepley 
1628514bf10dSMatthew G Knepley   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
16292d747510SLisandro Dalcin   if (pre && pre[0]) {
16302d747510SLisandro Dalcin     char *ptr = buf;
16319371c9d4SSatish Balay     if (name[0] == '-') {
16329371c9d4SSatish Balay       *ptr++ = '-';
16339371c9d4SSatish Balay       name++;
16349371c9d4SSatish Balay     }
16359b15cf9aSJacob Faibussowitsch     PetscCall(PetscStrncpy(ptr, pre, sizeof(buf) - ((ptr == buf) ? 0 : 1)));
16369566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
16372d747510SLisandro Dalcin     name = buf;
1638514bf10dSMatthew G Knepley   }
16392d747510SLisandro Dalcin 
164076bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
1641514bf10dSMatthew G Knepley     PetscBool valid;
16429355ec05SMatthew G. Knepley     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
16439566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
16449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(key, &valid));
164528b400f6SJacob Faibussowitsch     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1646514bf10dSMatthew G Knepley   }
1647514bf10dSMatthew G Knepley 
1648d6ced9c0SMatthew G. Knepley   /* determine the location and number of all _%d_ in the key */
1649d6ced9c0SMatthew G. Knepley   {
1650d6ced9c0SMatthew G. Knepley     int i, j;
1651d6ced9c0SMatthew G. Knepley     for (i = 0; name[i]; i++) {
1652d6ced9c0SMatthew G. Knepley       if (name[i] == '_') {
1653d6ced9c0SMatthew G. Knepley         for (j = i + 1; name[j]; j++) {
1654d6ced9c0SMatthew G. Knepley           if (name[j] >= '0' && name[j] <= '9') continue;
1655d6ced9c0SMatthew G. Knepley           if (name[j] == '_' && j > i + 1) { /* found a number */
1656d6ced9c0SMatthew G. Knepley             locs[numCnt]   = i + 1;
1657d6ced9c0SMatthew G. Knepley             loce[numCnt++] = j + 1;
1658d6ced9c0SMatthew G. Knepley           }
1659d6ced9c0SMatthew G. Knepley           i = j - 1;
1660d6ced9c0SMatthew G. Knepley           break;
1661d6ced9c0SMatthew G. Knepley         }
1662d6ced9c0SMatthew G. Knepley       }
1663d6ced9c0SMatthew G. Knepley     }
1664d6ced9c0SMatthew G. Knepley   }
1665d6ced9c0SMatthew G. Knepley 
1666363da2dcSJacob Faibussowitsch   /* slow search */
1667363da2dcSJacob Faibussowitsch   for (int c = -1; c < numCnt; ++c) {
1668363da2dcSJacob Faibussowitsch     char   opt[PETSC_MAX_OPTION_NAME + 2] = "";
16692d747510SLisandro Dalcin     size_t len;
1670d6ced9c0SMatthew G. Knepley 
1671d6ced9c0SMatthew G. Knepley     if (c < 0) {
1672c6a7a370SJeremy L Thompson       PetscCall(PetscStrncpy(opt, name, sizeof(opt)));
1673d6ced9c0SMatthew G. Knepley     } else {
1674363da2dcSJacob Faibussowitsch       PetscCall(PetscStrncpy(opt, name, PetscMin((size_t)(locs[c] + 1), sizeof(opt))));
1675363da2dcSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, name + loce[c], sizeof(opt) - 1));
1676d6ced9c0SMatthew G. Knepley     }
16779566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(opt, &len));
1678363da2dcSJacob Faibussowitsch     for (int i = 0; i < options->N; i++) {
1679363da2dcSJacob Faibussowitsch       PetscBool match;
1680363da2dcSJacob Faibussowitsch 
16819566063dSJacob Faibussowitsch       PetscCall(PetscStrncmp(options->names[i], opt, len, &match));
1682514bf10dSMatthew G Knepley       if (match) {
1683514bf10dSMatthew G Knepley         options->used[i] = PETSC_TRUE;
16842d747510SLisandro Dalcin         if (value) *value = options->values[i];
16852d747510SLisandro Dalcin         if (set) *set = PETSC_TRUE;
16863ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
1687514bf10dSMatthew G Knepley       }
1688514bf10dSMatthew G Knepley     }
16892d747510SLisandro Dalcin   }
16902d747510SLisandro Dalcin 
16912d747510SLisandro Dalcin   if (set) *set = PETSC_FALSE;
16923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1693514bf10dSMatthew G Knepley }
1694514bf10dSMatthew G Knepley 
1695e5c89e4eSSatish Balay /*@C
1696e5c89e4eSSatish Balay    PetscOptionsReject - Generates an error if a certain option is given.
1697e5c89e4eSSatish Balay 
16981c9f3c13SBarry Smith    Not Collective
1699e5c89e4eSSatish Balay 
1700e5c89e4eSSatish Balay    Input Parameters:
170120f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
170220f4b53cSBarry Smith .  pre - the option prefix (may be `NULL`)
17032d747510SLisandro Dalcin .  name - the option name one is seeking
170420f4b53cSBarry Smith -  mess - error message (may be `NULL`)
1705e5c89e4eSSatish Balay 
1706e5c89e4eSSatish Balay    Level: advanced
1707e5c89e4eSSatish Balay 
1708c2e3fba1SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1709db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1710db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1711c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1712db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1713db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1714e5c89e4eSSatish Balay @*/
1715d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsReject(PetscOptions options, const char pre[], const char name[], const char mess[])
1716d71ae5a4SJacob Faibussowitsch {
1717ace3abfcSBarry Smith   PetscBool flag = PETSC_FALSE;
1718e5c89e4eSSatish Balay 
1719e5c89e4eSSatish Balay   PetscFunctionBegin;
17209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(options, pre, name, &flag));
1721e5c89e4eSSatish Balay   if (flag) {
172208401ef6SPierre Jolivet     PetscCheck(!mess || !mess[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s with %s", pre ? pre : "", name + 1, mess);
1723f7d195e4SLawrence Mitchell     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s", pre ? pre : "", name + 1);
1724e5c89e4eSSatish Balay   }
17253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1726e5c89e4eSSatish Balay }
1727e5c89e4eSSatish Balay 
1728e5c89e4eSSatish Balay /*@C
17292d747510SLisandro Dalcin    PetscOptionsHasHelp - Determines whether the "-help" option is in the database.
17302d747510SLisandro Dalcin 
17312d747510SLisandro Dalcin    Not Collective
17322d747510SLisandro Dalcin 
1733811af0c4SBarry Smith    Input Parameter:
173420f4b53cSBarry Smith .  options - options database, use `NULL` for default global database
17352d747510SLisandro Dalcin 
1736811af0c4SBarry Smith    Output Parameter:
1737811af0c4SBarry Smith .  set - `PETSC_TRUE` if found else `PETSC_FALSE`.
17382d747510SLisandro Dalcin 
17392d747510SLisandro Dalcin    Level: advanced
17402d747510SLisandro Dalcin 
1741db781477SPatrick Sanan .seealso: `PetscOptionsHasName()`
17422d747510SLisandro Dalcin @*/
1743d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasHelp(PetscOptions options, PetscBool *set)
1744d71ae5a4SJacob Faibussowitsch {
17452d747510SLisandro Dalcin   PetscFunctionBegin;
1746dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(set, 2);
17472d747510SLisandro Dalcin   options = options ? options : defaultoptions;
17482d747510SLisandro Dalcin   *set    = options->help;
17493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17502d747510SLisandro Dalcin }
17512d747510SLisandro Dalcin 
1752d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options, PetscBool *set)
1753d71ae5a4SJacob Faibussowitsch {
1754d314f959SVaclav Hapla   PetscFunctionBegin;
1755dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(set, 2);
1756d314f959SVaclav Hapla   options = options ? options : defaultoptions;
1757d314f959SVaclav Hapla   *set    = options->help_intro;
17583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1759d314f959SVaclav Hapla }
1760d314f959SVaclav Hapla 
17612d747510SLisandro Dalcin /*@C
1762e24fcbf7SPierre Jolivet    PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or Boolean, even
1763e24fcbf7SPierre Jolivet                       if its value is set to false.
1764e5c89e4eSSatish Balay 
1765e5c89e4eSSatish Balay    Not Collective
1766e5c89e4eSSatish Balay 
1767e5c89e4eSSatish Balay    Input Parameters:
176820f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
176920f4b53cSBarry Smith .  pre - string to prepend to the name or `NULL`
17703de71b31SHong Zhang -  name - the option one is seeking
1771e5c89e4eSSatish Balay 
1772811af0c4SBarry Smith    Output Parameter:
1773811af0c4SBarry Smith .  set - `PETSC_TRUE` if found else `PETSC_FALSE`.
1774e5c89e4eSSatish Balay 
1775e5c89e4eSSatish Balay    Level: beginner
1776e5c89e4eSSatish Balay 
1777811af0c4SBarry Smith    Note:
1778811af0c4SBarry Smith    In many cases you probably want to use `PetscOptionsGetBool()` instead of calling this, to allowing toggling values.
177990d69ab7SBarry Smith 
1780db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
1781db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1782db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1783c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1784db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1785db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1786e5c89e4eSSatish Balay @*/
1787d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasName(PetscOptions options, const char pre[], const char name[], PetscBool *set)
1788d71ae5a4SJacob Faibussowitsch {
17892d747510SLisandro Dalcin   const char *value;
1790ace3abfcSBarry Smith   PetscBool   flag;
1791e5c89e4eSSatish Balay 
1792e5c89e4eSSatish Balay   PetscFunctionBegin;
17939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
179496ef3cdfSSatish Balay   if (set) *set = flag;
17953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1796e5c89e4eSSatish Balay }
1797e5c89e4eSSatish Balay 
1798e5c89e4eSSatish Balay /*@C
17992d747510SLisandro Dalcin    PetscOptionsGetAll - Lists all the options the program was run with in a single string.
18002d747510SLisandro Dalcin 
18012d747510SLisandro Dalcin    Not Collective
18022d747510SLisandro Dalcin 
1803fd292e60Sprj-    Input Parameter:
180420f4b53cSBarry Smith .  options - the options database, use `NULL` for the default global database
18052d747510SLisandro Dalcin 
18062d747510SLisandro Dalcin    Output Parameter:
18072d747510SLisandro Dalcin .  copts - pointer where string pointer is stored
18082d747510SLisandro Dalcin 
180920f4b53cSBarry Smith    Level: advanced
181020f4b53cSBarry Smith 
18112d747510SLisandro Dalcin    Notes:
1812811af0c4SBarry Smith     The array and each entry in the array should be freed with `PetscFree()`
1813811af0c4SBarry Smith 
18141c9f3c13SBarry Smith     Each process may have different values depending on how the options were inserted into the database
18152d747510SLisandro Dalcin 
1816db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsView()`, `PetscOptionsPush()`, `PetscOptionsPop()`
18172d747510SLisandro Dalcin @*/
1818d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetAll(PetscOptions options, char *copts[])
1819d71ae5a4SJacob Faibussowitsch {
18202d747510SLisandro Dalcin   PetscInt i;
18212d747510SLisandro Dalcin   size_t   len = 1, lent = 0;
18222d747510SLisandro Dalcin   char    *coptions = NULL;
18232d747510SLisandro Dalcin 
18242d747510SLisandro Dalcin   PetscFunctionBegin;
18252d747510SLisandro Dalcin   PetscValidPointer(copts, 2);
18262d747510SLisandro Dalcin   options = options ? options : defaultoptions;
18272d747510SLisandro Dalcin   /* count the length of the required string */
18282d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
18299566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(options->names[i], &lent));
18302d747510SLisandro Dalcin     len += 2 + lent;
18312d747510SLisandro Dalcin     if (options->values[i]) {
18329566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(options->values[i], &lent));
18332d747510SLisandro Dalcin       len += 1 + lent;
18342d747510SLisandro Dalcin     }
18352d747510SLisandro Dalcin   }
18369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &coptions));
18372d747510SLisandro Dalcin   coptions[0] = 0;
18382d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
1839c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, "-", len));
1840c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, options->names[i], len));
1841c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, " ", len));
18422d747510SLisandro Dalcin     if (options->values[i]) {
1843c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(coptions, options->values[i], len));
1844c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(coptions, " ", len));
18452d747510SLisandro Dalcin     }
18462d747510SLisandro Dalcin   }
18472d747510SLisandro Dalcin   *copts = coptions;
18483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18492d747510SLisandro Dalcin }
18502d747510SLisandro Dalcin 
18512d747510SLisandro Dalcin /*@C
18522d747510SLisandro Dalcin    PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database
18532d747510SLisandro Dalcin 
18542d747510SLisandro Dalcin    Not Collective
18552d747510SLisandro Dalcin 
1856d8d19677SJose E. Roman    Input Parameters:
185720f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
18582d747510SLisandro Dalcin -  name - string name of option
18592d747510SLisandro Dalcin 
18602d747510SLisandro Dalcin    Output Parameter:
1861811af0c4SBarry Smith .  used - `PETSC_TRUE` if the option was used, otherwise false, including if option was not found in options database
18622d747510SLisandro Dalcin 
18632d747510SLisandro Dalcin    Level: advanced
18642d747510SLisandro Dalcin 
1865811af0c4SBarry Smith    Note:
18669666a313SBarry Smith    The value returned may be different on each process and depends on which options have been processed
18671c9f3c13SBarry Smith    on the given process
18681c9f3c13SBarry Smith 
1869db781477SPatrick Sanan .seealso: `PetscOptionsView()`, `PetscOptionsLeft()`, `PetscOptionsAllUsed()`
18702d747510SLisandro Dalcin @*/
1871d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsUsed(PetscOptions options, const char *name, PetscBool *used)
1872d71ae5a4SJacob Faibussowitsch {
18732d747510SLisandro Dalcin   PetscInt i;
18742d747510SLisandro Dalcin 
18752d747510SLisandro Dalcin   PetscFunctionBegin;
18762d747510SLisandro Dalcin   PetscValidCharPointer(name, 2);
1877dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(used, 3);
18782d747510SLisandro Dalcin   options = options ? options : defaultoptions;
18792d747510SLisandro Dalcin   *used   = PETSC_FALSE;
18802d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
18819566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(options->names[i], name, used));
18822d747510SLisandro Dalcin     if (*used) {
18832d747510SLisandro Dalcin       *used = options->used[i];
18842d747510SLisandro Dalcin       break;
18852d747510SLisandro Dalcin     }
18862d747510SLisandro Dalcin   }
18873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18882d747510SLisandro Dalcin }
18892d747510SLisandro Dalcin 
1890487a658cSBarry Smith /*@
18912d747510SLisandro Dalcin    PetscOptionsAllUsed - Returns a count of the number of options in the
18922d747510SLisandro Dalcin    database that have never been selected.
18932d747510SLisandro Dalcin 
18942d747510SLisandro Dalcin    Not Collective
18952d747510SLisandro Dalcin 
18962d747510SLisandro Dalcin    Input Parameter:
189720f4b53cSBarry Smith .  options - options database, use `NULL` for default global database
18982d747510SLisandro Dalcin 
18992d747510SLisandro Dalcin    Output Parameter:
19002d747510SLisandro Dalcin .  N - count of options not used
19012d747510SLisandro Dalcin 
19022d747510SLisandro Dalcin    Level: advanced
19032d747510SLisandro Dalcin 
1904811af0c4SBarry Smith    Note:
19059666a313SBarry Smith    The value returned may be different on each process and depends on which options have been processed
19061c9f3c13SBarry Smith    on the given process
19071c9f3c13SBarry Smith 
1908db781477SPatrick Sanan .seealso: `PetscOptionsView()`
19092d747510SLisandro Dalcin @*/
1910d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsAllUsed(PetscOptions options, PetscInt *N)
1911d71ae5a4SJacob Faibussowitsch {
19122d747510SLisandro Dalcin   PetscInt i, n = 0;
19132d747510SLisandro Dalcin 
19142d747510SLisandro Dalcin   PetscFunctionBegin;
19152d747510SLisandro Dalcin   PetscValidIntPointer(N, 2);
19162d747510SLisandro Dalcin   options = options ? options : defaultoptions;
19172d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
19182d747510SLisandro Dalcin     if (!options->used[i]) n++;
19192d747510SLisandro Dalcin   }
19202d747510SLisandro Dalcin   *N = n;
19213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19222d747510SLisandro Dalcin }
19232d747510SLisandro Dalcin 
1924487a658cSBarry Smith /*@
19252d747510SLisandro Dalcin    PetscOptionsLeft - Prints to screen any options that were set and never used.
19262d747510SLisandro Dalcin 
19272d747510SLisandro Dalcin    Not Collective
19282d747510SLisandro Dalcin 
19292d747510SLisandro Dalcin    Input Parameter:
193020f4b53cSBarry Smith .  options - options database; use `NULL` for default global database
19312d747510SLisandro Dalcin 
19322d747510SLisandro Dalcin    Options Database Key:
1933811af0c4SBarry Smith .  -options_left - activates `PetscOptionsAllUsed()` within `PetscFinalize()`
19342d747510SLisandro Dalcin 
193520f4b53cSBarry Smith    Level: advanced
193620f4b53cSBarry Smith 
19373de2bfdfSBarry Smith    Notes:
1938811af0c4SBarry Smith       This is rarely used directly, it is called by `PetscFinalize()` in debug more or if -options_left
19391c9f3c13SBarry Smith       is passed otherwise to help users determine possible mistakes in their usage of options. This
1940811af0c4SBarry Smith       only prints values on process zero of `PETSC_COMM_WORLD`.
1941811af0c4SBarry Smith 
1942811af0c4SBarry Smith       Other processes depending the objects
19431c9f3c13SBarry Smith       used may have different options that are left unused.
19443de2bfdfSBarry Smith 
1945db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`
19462d747510SLisandro Dalcin @*/
1947d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeft(PetscOptions options)
1948d71ae5a4SJacob Faibussowitsch {
19492d747510SLisandro Dalcin   PetscInt     i;
19503de2bfdfSBarry Smith   PetscInt     cnt = 0;
19513de2bfdfSBarry Smith   PetscOptions toptions;
19522d747510SLisandro Dalcin 
19532d747510SLisandro Dalcin   PetscFunctionBegin;
19543de2bfdfSBarry Smith   toptions = options ? options : defaultoptions;
19553de2bfdfSBarry Smith   for (i = 0; i < toptions->N; i++) {
19563de2bfdfSBarry Smith     if (!toptions->used[i]) {
1957660278c0SBarry Smith       if (PetscCIOption(toptions->names[i])) continue;
19583de2bfdfSBarry Smith       if (toptions->values[i]) {
19599355ec05SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s value: %s source: %s\n", toptions->names[i], toptions->values[i], PetscOptionSources[toptions->source[i]]));
19602d747510SLisandro Dalcin       } else {
19619355ec05SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s (no value) source: %s\n", toptions->names[i], PetscOptionSources[toptions->source[i]]));
19622d747510SLisandro Dalcin       }
19632d747510SLisandro Dalcin     }
19642d747510SLisandro Dalcin   }
19653de2bfdfSBarry Smith   if (!options) {
19663de2bfdfSBarry Smith     toptions = defaultoptions;
19673de2bfdfSBarry Smith     while (toptions->previous) {
19683de2bfdfSBarry Smith       cnt++;
19693de2bfdfSBarry Smith       toptions = toptions->previous;
19703de2bfdfSBarry Smith     }
197148a46eb9SPierre Jolivet     if (cnt) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: You may have forgotten some calls to PetscOptionsPop(),\n             PetscOptionsPop() has been called %" PetscInt_FMT " less times than PetscOptionsPush()\n", cnt));
19723de2bfdfSBarry Smith   }
19733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19742d747510SLisandro Dalcin }
19752d747510SLisandro Dalcin 
19762d747510SLisandro Dalcin /*@C
19772d747510SLisandro Dalcin    PetscOptionsLeftGet - Returns all options that were set and never used.
19782d747510SLisandro Dalcin 
19792d747510SLisandro Dalcin    Not Collective
19802d747510SLisandro Dalcin 
19812d747510SLisandro Dalcin    Input Parameter:
198220f4b53cSBarry Smith .  options - options database, use `NULL` for default global database
19832d747510SLisandro Dalcin 
1984d8d19677SJose E. Roman    Output Parameters:
1985a2b725a8SWilliam Gropp +  N - count of options not used
19862d747510SLisandro Dalcin .  names - names of options not used
1987a2b725a8SWilliam Gropp -  values - values of options not used
19882d747510SLisandro Dalcin 
19892d747510SLisandro Dalcin    Level: advanced
19902d747510SLisandro Dalcin 
19912d747510SLisandro Dalcin    Notes:
1992811af0c4SBarry Smith    Users should call `PetscOptionsLeftRestore()` to free the memory allocated in this routine
1993811af0c4SBarry Smith 
1994811af0c4SBarry Smith    The value returned may be different on each process and depends on which options have been processed
19951c9f3c13SBarry Smith    on the given process
19962d747510SLisandro Dalcin 
1997db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`
19982d747510SLisandro Dalcin @*/
1999d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftGet(PetscOptions options, PetscInt *N, char **names[], char **values[])
2000d71ae5a4SJacob Faibussowitsch {
20012d747510SLisandro Dalcin   PetscInt i, n;
20022d747510SLisandro Dalcin 
20032d747510SLisandro Dalcin   PetscFunctionBegin;
20042d747510SLisandro Dalcin   if (N) PetscValidIntPointer(N, 2);
20052d747510SLisandro Dalcin   if (names) PetscValidPointer(names, 3);
20062d747510SLisandro Dalcin   if (values) PetscValidPointer(values, 4);
20072d747510SLisandro Dalcin   options = options ? options : defaultoptions;
20082d747510SLisandro Dalcin 
20092d747510SLisandro Dalcin   /* The number of unused PETSc options */
20102d747510SLisandro Dalcin   n = 0;
20112d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
2012660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
20132d747510SLisandro Dalcin     if (!options->used[i]) n++;
20142d747510SLisandro Dalcin   }
2015ad540459SPierre Jolivet   if (N) *N = n;
20169566063dSJacob Faibussowitsch   if (names) PetscCall(PetscMalloc1(n, names));
20179566063dSJacob Faibussowitsch   if (values) PetscCall(PetscMalloc1(n, values));
20182d747510SLisandro Dalcin 
20192d747510SLisandro Dalcin   n = 0;
20202d747510SLisandro Dalcin   if (names || values) {
20212d747510SLisandro Dalcin     for (i = 0; i < options->N; i++) {
20222d747510SLisandro Dalcin       if (!options->used[i]) {
2023660278c0SBarry Smith         if (PetscCIOption(options->names[i])) continue;
20242d747510SLisandro Dalcin         if (names) (*names)[n] = options->names[i];
20252d747510SLisandro Dalcin         if (values) (*values)[n] = options->values[i];
20262d747510SLisandro Dalcin         n++;
20272d747510SLisandro Dalcin       }
20282d747510SLisandro Dalcin     }
20292d747510SLisandro Dalcin   }
20303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20312d747510SLisandro Dalcin }
20322d747510SLisandro Dalcin 
20332d747510SLisandro Dalcin /*@C
2034811af0c4SBarry Smith    PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using `PetscOptionsLeftGet()`.
20352d747510SLisandro Dalcin 
20362d747510SLisandro Dalcin    Not Collective
20372d747510SLisandro Dalcin 
2038d8d19677SJose E. Roman    Input Parameters:
203920f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
20402d747510SLisandro Dalcin .  names - names of options not used
2041a2b725a8SWilliam Gropp -  values - values of options not used
20422d747510SLisandro Dalcin 
20432d747510SLisandro Dalcin    Level: advanced
20442d747510SLisandro Dalcin 
2045db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`, `PetscOptionsLeftGet()`
20462d747510SLisandro Dalcin @*/
2047d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftRestore(PetscOptions options, PetscInt *N, char **names[], char **values[])
2048d71ae5a4SJacob Faibussowitsch {
20492d747510SLisandro Dalcin   PetscFunctionBegin;
20502d747510SLisandro Dalcin   if (N) PetscValidIntPointer(N, 2);
20512d747510SLisandro Dalcin   if (names) PetscValidPointer(names, 3);
20522d747510SLisandro Dalcin   if (values) PetscValidPointer(values, 4);
2053ad540459SPierre Jolivet   if (N) *N = 0;
20549566063dSJacob Faibussowitsch   if (names) PetscCall(PetscFree(*names));
20559566063dSJacob Faibussowitsch   if (values) PetscCall(PetscFree(*values));
20563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20572d747510SLisandro Dalcin }
20582d747510SLisandro Dalcin 
20592d747510SLisandro Dalcin /*@C
2060811af0c4SBarry Smith    PetscOptionsMonitorDefault - Print all options set value events using the supplied `PetscViewer`.
20612d747510SLisandro Dalcin 
2062c3339decSBarry Smith    Logically Collective
20632d747510SLisandro Dalcin 
20642d747510SLisandro Dalcin    Input Parameters:
20652d747510SLisandro Dalcin +  name  - option name string
20662d747510SLisandro Dalcin .  value - option value string
20679355ec05SMatthew G. Knepley .  source - The source for the option
206820f4b53cSBarry Smith -  ctx - a `PETSCVIEWERASCII` or `NULL`
20692d747510SLisandro Dalcin 
20702d747510SLisandro Dalcin    Level: intermediate
20712d747510SLisandro Dalcin 
20729666a313SBarry Smith    Notes:
207320f4b53cSBarry Smith      If ctx is `NULL`, `PetscPrintf()` is used.
2074811af0c4SBarry Smith      The first MPI rank in the `PetscViewer` viewer actually prints the values, other
20751c9f3c13SBarry Smith      processes may have different values set
20761c9f3c13SBarry Smith 
2077811af0c4SBarry Smith      If `PetscCIEnabled` then do not print the test harness options
2078660278c0SBarry Smith 
2079db781477SPatrick Sanan .seealso: `PetscOptionsMonitorSet()`
20802d747510SLisandro Dalcin @*/
20819355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsMonitorDefault(const char name[], const char value[], PetscOptionSource source, void *ctx)
2082d71ae5a4SJacob Faibussowitsch {
20832d747510SLisandro Dalcin   PetscFunctionBegin;
20843ba16761SJacob Faibussowitsch   if (PetscCIOption(name)) PetscFunctionReturn(PETSC_SUCCESS);
2085660278c0SBarry Smith 
20869060e2f9SVaclav Hapla   if (ctx) {
20879060e2f9SVaclav Hapla     PetscViewer viewer = (PetscViewer)ctx;
20882d747510SLisandro Dalcin     if (!value) {
20899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Removing option: %s\n", name));
20902d747510SLisandro Dalcin     } else if (!value[0]) {
20919355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
20922d747510SLisandro Dalcin     } else {
20939355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
20942d747510SLisandro Dalcin     }
20959060e2f9SVaclav Hapla   } else {
20969060e2f9SVaclav Hapla     MPI_Comm comm = PETSC_COMM_WORLD;
20979060e2f9SVaclav Hapla     if (!value) {
20989566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Removing option: %s\n", name));
20999060e2f9SVaclav Hapla     } else if (!value[0]) {
21009355ec05SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
21019060e2f9SVaclav Hapla     } else {
2102aaa8cc7dSPierre Jolivet       PetscCall(PetscPrintf(comm, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
21039060e2f9SVaclav Hapla     }
21049060e2f9SVaclav Hapla   }
21053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21062d747510SLisandro Dalcin }
21072d747510SLisandro Dalcin 
21082d747510SLisandro Dalcin /*@C
21092d747510SLisandro Dalcin    PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
21102d747510SLisandro Dalcin    modified the PETSc options database.
21112d747510SLisandro Dalcin 
21122d747510SLisandro Dalcin    Not Collective
21132d747510SLisandro Dalcin 
21142d747510SLisandro Dalcin    Input Parameters:
211520f4b53cSBarry Smith +  monitor - pointer to function (if this is `NULL`, it turns off monitoring
21162d747510SLisandro Dalcin .  mctx    - [optional] context for private data for the
211720f4b53cSBarry Smith              monitor routine (use `NULL` if no context is desired)
21182d747510SLisandro Dalcin -  monitordestroy - [optional] routine that frees monitor context
211920f4b53cSBarry Smith           (may be `NULL`)
21202d747510SLisandro Dalcin 
212120f4b53cSBarry Smith    Calling Sequence of `monitor`:
212220f4b53cSBarry Smith $   PetscErrorCode monitor(const char name[], const char value[], void *mctx)
21232d747510SLisandro Dalcin +  name - option name string
21242d747510SLisandro Dalcin .  value - option value string
21259355ec05SMatthew G. Knepley . source - option source
2126811af0c4SBarry Smith -  mctx  - optional monitoring context, as set by `PetscOptionsMonitorSet()`
21272d747510SLisandro Dalcin 
212820f4b53cSBarry Smith    Calling Sequence of `monitordestroy`:
212920f4b53cSBarry Smith $  PetscErrorCode monitordestroy(void *cctx)
213020f4b53cSBarry Smith 
213120f4b53cSBarry Smith    Options Database Key:
2132811af0c4SBarry Smith    See `PetscInitialize()` for options related to option database monitoring.
21332d747510SLisandro Dalcin 
213420f4b53cSBarry Smith    Level: intermediate
213520f4b53cSBarry Smith 
21362d747510SLisandro Dalcin    Notes:
21372d747510SLisandro Dalcin    The default is to do nothing.  To print the name and value of options
2138811af0c4SBarry Smith    being inserted into the database, use `PetscOptionsMonitorDefault()` as the monitoring routine,
21392d747510SLisandro Dalcin    with a null monitoring context.
21402d747510SLisandro Dalcin 
21412d747510SLisandro Dalcin    Several different monitoring routines may be set by calling
2142811af0c4SBarry Smith    `PetscOptionsMonitorSet()` multiple times; all will be called in the
21432d747510SLisandro Dalcin    order in which they were set.
21442d747510SLisandro Dalcin 
2145db781477SPatrick Sanan .seealso: `PetscOptionsMonitorDefault()`, `PetscInitialize()`
21462d747510SLisandro Dalcin @*/
21479355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], PetscOptionSource, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
2148d71ae5a4SJacob Faibussowitsch {
21492d747510SLisandro Dalcin   PetscOptions options = defaultoptions;
21502d747510SLisandro Dalcin 
21512d747510SLisandro Dalcin   PetscFunctionBegin;
21523ba16761SJacob Faibussowitsch   if (options->monitorCancel) PetscFunctionReturn(PETSC_SUCCESS);
215308401ef6SPierre Jolivet   PetscCheck(options->numbermonitors < MAXOPTIONSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many PetscOptions monitors set");
21542d747510SLisandro Dalcin   options->monitor[options->numbermonitors]          = monitor;
21552d747510SLisandro Dalcin   options->monitordestroy[options->numbermonitors]   = monitordestroy;
21562d747510SLisandro Dalcin   options->monitorcontext[options->numbermonitors++] = (void *)mctx;
21573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21582d747510SLisandro Dalcin }
21592d747510SLisandro Dalcin 
21602d747510SLisandro Dalcin /*
21612d747510SLisandro Dalcin    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
216263fe8743SVaclav Hapla      Empty string is considered as true.
21632d747510SLisandro Dalcin */
2164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToBool(const char value[], PetscBool *a)
2165d71ae5a4SJacob Faibussowitsch {
21662d747510SLisandro Dalcin   PetscBool istrue, isfalse;
21672d747510SLisandro Dalcin   size_t    len;
21682d747510SLisandro Dalcin 
21692d747510SLisandro Dalcin   PetscFunctionBegin;
217063fe8743SVaclav Hapla   /* PetscStrlen() returns 0 for NULL or "" */
21719566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(value, &len));
21729371c9d4SSatish Balay   if (!len) {
21739371c9d4SSatish Balay     *a = PETSC_TRUE;
21743ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21759371c9d4SSatish Balay   }
21769566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "TRUE", &istrue));
21779371c9d4SSatish Balay   if (istrue) {
21789371c9d4SSatish Balay     *a = PETSC_TRUE;
21793ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21809371c9d4SSatish Balay   }
21819566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "YES", &istrue));
21829371c9d4SSatish Balay   if (istrue) {
21839371c9d4SSatish Balay     *a = PETSC_TRUE;
21843ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21859371c9d4SSatish Balay   }
21869566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "1", &istrue));
21879371c9d4SSatish Balay   if (istrue) {
21889371c9d4SSatish Balay     *a = PETSC_TRUE;
21893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21909371c9d4SSatish Balay   }
21919566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "on", &istrue));
21929371c9d4SSatish Balay   if (istrue) {
21939371c9d4SSatish Balay     *a = PETSC_TRUE;
21943ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21959371c9d4SSatish Balay   }
21969566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "FALSE", &isfalse));
21979371c9d4SSatish Balay   if (isfalse) {
21989371c9d4SSatish Balay     *a = PETSC_FALSE;
21993ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22009371c9d4SSatish Balay   }
22019566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "NO", &isfalse));
22029371c9d4SSatish Balay   if (isfalse) {
22039371c9d4SSatish Balay     *a = PETSC_FALSE;
22043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22059371c9d4SSatish Balay   }
22069566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "0", &isfalse));
22079371c9d4SSatish Balay   if (isfalse) {
22089371c9d4SSatish Balay     *a = PETSC_FALSE;
22093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22109371c9d4SSatish Balay   }
22119566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "off", &isfalse));
22129371c9d4SSatish Balay   if (isfalse) {
22139371c9d4SSatish Balay     *a = PETSC_FALSE;
22143ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22159371c9d4SSatish Balay   }
221698921bdaSJacob Faibussowitsch   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown logical value: %s", value);
22172d747510SLisandro Dalcin }
22182d747510SLisandro Dalcin 
22192d747510SLisandro Dalcin /*
22202d747510SLisandro Dalcin    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
22212d747510SLisandro Dalcin */
2222d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToInt(const char name[], PetscInt *a)
2223d71ae5a4SJacob Faibussowitsch {
22242d747510SLisandro Dalcin   size_t    len;
22252d747510SLisandro Dalcin   PetscBool decide, tdefault, mouse;
22262d747510SLisandro Dalcin 
22272d747510SLisandro Dalcin   PetscFunctionBegin;
22289566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
22295f80ce2aSJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
22302d747510SLisandro Dalcin 
22319566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &tdefault));
223248a46eb9SPierre Jolivet   if (!tdefault) PetscCall(PetscStrcasecmp(name, "DEFAULT", &tdefault));
22339566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &decide));
223448a46eb9SPierre Jolivet   if (!decide) PetscCall(PetscStrcasecmp(name, "DECIDE", &decide));
22359566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "mouse", &mouse));
22362d747510SLisandro Dalcin 
22372d747510SLisandro Dalcin   if (tdefault) *a = PETSC_DEFAULT;
22382d747510SLisandro Dalcin   else if (decide) *a = PETSC_DECIDE;
22392d747510SLisandro Dalcin   else if (mouse) *a = -1;
22402d747510SLisandro Dalcin   else {
22412d747510SLisandro Dalcin     char *endptr;
22422d747510SLisandro Dalcin     long  strtolval;
22432d747510SLisandro Dalcin 
22442d747510SLisandro Dalcin     strtolval = strtol(name, &endptr, 10);
2245cc73adaaSBarry Smith     PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no integer value (do not include . in it)", name);
22462d747510SLisandro Dalcin 
22472d747510SLisandro Dalcin #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
22482d747510SLisandro Dalcin     (void)strtolval;
22492d747510SLisandro Dalcin     *a = atoll(name);
22502d747510SLisandro Dalcin #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
22512d747510SLisandro Dalcin     (void)strtolval;
22522d747510SLisandro Dalcin     *a = _atoi64(name);
22532d747510SLisandro Dalcin #else
22542d747510SLisandro Dalcin     *a = (PetscInt)strtolval;
22552d747510SLisandro Dalcin #endif
22562d747510SLisandro Dalcin   }
22573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22582d747510SLisandro Dalcin }
22592d747510SLisandro Dalcin 
22602d747510SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
22612d747510SLisandro Dalcin   #include <quadmath.h>
22622d747510SLisandro Dalcin #endif
22632d747510SLisandro Dalcin 
2264d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscStrtod(const char name[], PetscReal *a, char **endptr)
2265d71ae5a4SJacob Faibussowitsch {
22662d747510SLisandro Dalcin   PetscFunctionBegin;
22672d747510SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
22682d747510SLisandro Dalcin   *a = strtoflt128(name, endptr);
22692d747510SLisandro Dalcin #else
22702d747510SLisandro Dalcin   *a = (PetscReal)strtod(name, endptr);
22712d747510SLisandro Dalcin #endif
22723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22732d747510SLisandro Dalcin }
22742d747510SLisandro Dalcin 
2275d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscStrtoz(const char name[], PetscScalar *a, char **endptr, PetscBool *isImaginary)
2276d71ae5a4SJacob Faibussowitsch {
22772d747510SLisandro Dalcin   PetscBool hasi = PETSC_FALSE;
22782d747510SLisandro Dalcin   char     *ptr;
22792d747510SLisandro Dalcin   PetscReal strtoval;
22802d747510SLisandro Dalcin 
22812d747510SLisandro Dalcin   PetscFunctionBegin;
22829566063dSJacob Faibussowitsch   PetscCall(PetscStrtod(name, &strtoval, &ptr));
22832d747510SLisandro Dalcin   if (ptr == name) {
22842d747510SLisandro Dalcin     strtoval = 1.;
22852d747510SLisandro Dalcin     hasi     = PETSC_TRUE;
22862d747510SLisandro Dalcin     if (name[0] == 'i') {
22872d747510SLisandro Dalcin       ptr++;
22882d747510SLisandro Dalcin     } else if (name[0] == '+' && name[1] == 'i') {
22892d747510SLisandro Dalcin       ptr += 2;
22902d747510SLisandro Dalcin     } else if (name[0] == '-' && name[1] == 'i') {
22912d747510SLisandro Dalcin       strtoval = -1.;
22922d747510SLisandro Dalcin       ptr += 2;
22932d747510SLisandro Dalcin     }
22942d747510SLisandro Dalcin   } else if (*ptr == 'i') {
22952d747510SLisandro Dalcin     hasi = PETSC_TRUE;
22962d747510SLisandro Dalcin     ptr++;
22972d747510SLisandro Dalcin   }
22982d747510SLisandro Dalcin   *endptr      = ptr;
22992d747510SLisandro Dalcin   *isImaginary = hasi;
23002d747510SLisandro Dalcin   if (hasi) {
23012d747510SLisandro Dalcin #if !defined(PETSC_USE_COMPLEX)
230298921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s contains imaginary but complex not supported ", name);
23032d747510SLisandro Dalcin #else
23042d747510SLisandro Dalcin     *a = PetscCMPLX(0., strtoval);
23052d747510SLisandro Dalcin #endif
23062d747510SLisandro Dalcin   } else {
23072d747510SLisandro Dalcin     *a = strtoval;
23082d747510SLisandro Dalcin   }
23093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23102d747510SLisandro Dalcin }
23112d747510SLisandro Dalcin 
23122d747510SLisandro Dalcin /*
23132d747510SLisandro Dalcin    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
23142d747510SLisandro Dalcin */
2315d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToReal(const char name[], PetscReal *a)
2316d71ae5a4SJacob Faibussowitsch {
23172d747510SLisandro Dalcin   size_t    len;
23182d747510SLisandro Dalcin   PetscBool match;
23192d747510SLisandro Dalcin   char     *endptr;
23202d747510SLisandro Dalcin 
23212d747510SLisandro Dalcin   PetscFunctionBegin;
23229566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
232328b400f6SJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "String of length zero has no numerical value");
23242d747510SLisandro Dalcin 
23259566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &match));
23269566063dSJacob Faibussowitsch   if (!match) PetscCall(PetscStrcasecmp(name, "DEFAULT", &match));
23279371c9d4SSatish Balay   if (match) {
23289371c9d4SSatish Balay     *a = PETSC_DEFAULT;
23293ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
23309371c9d4SSatish Balay   }
23312d747510SLisandro Dalcin 
23329566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &match));
23339566063dSJacob Faibussowitsch   if (!match) PetscCall(PetscStrcasecmp(name, "DECIDE", &match));
23349371c9d4SSatish Balay   if (match) {
23359371c9d4SSatish Balay     *a = PETSC_DECIDE;
23363ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
23379371c9d4SSatish Balay   }
23382d747510SLisandro Dalcin 
23399566063dSJacob Faibussowitsch   PetscCall(PetscStrtod(name, a, &endptr));
234039a651e2SJacob Faibussowitsch   PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value", name);
23413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23422d747510SLisandro Dalcin }
23432d747510SLisandro Dalcin 
2344d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToScalar(const char name[], PetscScalar *a)
2345d71ae5a4SJacob Faibussowitsch {
23462d747510SLisandro Dalcin   PetscBool   imag1;
23472d747510SLisandro Dalcin   size_t      len;
23482d747510SLisandro Dalcin   PetscScalar val = 0.;
23492d747510SLisandro Dalcin   char       *ptr = NULL;
23502d747510SLisandro Dalcin 
23512d747510SLisandro Dalcin   PetscFunctionBegin;
23529566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
235328b400f6SJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
23549566063dSJacob Faibussowitsch   PetscCall(PetscStrtoz(name, &val, &ptr, &imag1));
23552d747510SLisandro Dalcin #if defined(PETSC_USE_COMPLEX)
23562d747510SLisandro Dalcin   if ((size_t)(ptr - name) < len) {
23572d747510SLisandro Dalcin     PetscBool   imag2;
23582d747510SLisandro Dalcin     PetscScalar val2;
23592d747510SLisandro Dalcin 
23609566063dSJacob Faibussowitsch     PetscCall(PetscStrtoz(ptr, &val2, &ptr, &imag2));
236139a651e2SJacob Faibussowitsch     if (imag1) PetscCheck(imag2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s: must specify imaginary component second", name);
23622d747510SLisandro Dalcin     val = PetscCMPLX(PetscRealPart(val), PetscImaginaryPart(val2));
23632d747510SLisandro Dalcin   }
23642d747510SLisandro Dalcin #endif
236539a651e2SJacob Faibussowitsch   PetscCheck((size_t)(ptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value ", name);
23662d747510SLisandro Dalcin   *a = val;
23673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23682d747510SLisandro Dalcin }
23692d747510SLisandro Dalcin 
23702d747510SLisandro Dalcin /*@C
23712d747510SLisandro Dalcin    PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
23722d747510SLisandro Dalcin             option in the database.
2373e5c89e4eSSatish Balay 
2374e5c89e4eSSatish Balay    Not Collective
2375e5c89e4eSSatish Balay 
2376e5c89e4eSSatish Balay    Input Parameters:
237720f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
237820f4b53cSBarry Smith .  pre - the string to prepend to the name or `NULL`
2379e5c89e4eSSatish Balay -  name - the option one is seeking
2380e5c89e4eSSatish Balay 
2381d8d19677SJose E. Roman    Output Parameters:
23822d747510SLisandro Dalcin +  ivalue - the logical value to return
2383811af0c4SBarry Smith -  set - `PETSC_TRUE`  if found, else `PETSC_FALSE`
2384e5c89e4eSSatish Balay 
2385e5c89e4eSSatish Balay    Level: beginner
2386e5c89e4eSSatish Balay 
238795452b02SPatrick Sanan    Notes:
2388811af0c4SBarry Smith        TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
2389811af0c4SBarry Smith        FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
23902d747510SLisandro Dalcin 
2391811af0c4SBarry Smith       If the option is given, but no value is provided, then ivalue and set are both given the value `PETSC_TRUE`. That is -requested_bool
23922d747510SLisandro Dalcin      is equivalent to -requested_bool true
23932d747510SLisandro Dalcin 
23942d747510SLisandro Dalcin        If the user does not supply the option at all ivalue is NOT changed. Thus
23952efd9cb1SBarry Smith      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
23962efd9cb1SBarry Smith 
2397db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2398db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2399db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2400c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2401db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2402db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2403e5c89e4eSSatish Balay @*/
2404d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetBool(PetscOptions options, const char pre[], const char name[], PetscBool *ivalue, PetscBool *set)
2405d71ae5a4SJacob Faibussowitsch {
24062d747510SLisandro Dalcin   const char *value;
2407ace3abfcSBarry Smith   PetscBool   flag;
2408e5c89e4eSSatish Balay 
2409e5c89e4eSSatish Balay   PetscFunctionBegin;
24102d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
2411064a246eSJacob Faibussowitsch   if (ivalue) PetscValidBoolPointer(ivalue, 4);
24129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2413e5c89e4eSSatish Balay   if (flag) {
241496ef3cdfSSatish Balay     if (set) *set = PETSC_TRUE;
24159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(value, &flag));
24162d747510SLisandro Dalcin     if (ivalue) *ivalue = flag;
2417e5c89e4eSSatish Balay   } else {
241896ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2419e5c89e4eSSatish Balay   }
24203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2421e5c89e4eSSatish Balay }
2422e5c89e4eSSatish Balay 
2423e5c89e4eSSatish Balay /*@C
2424e5c89e4eSSatish Balay    PetscOptionsGetEList - Puts a list of option values that a single one may be selected from
2425e5c89e4eSSatish Balay 
2426e5c89e4eSSatish Balay    Not Collective
2427e5c89e4eSSatish Balay 
2428e5c89e4eSSatish Balay    Input Parameters:
242920f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
243020f4b53cSBarry Smith .  pre - the string to prepend to the name or `NULL`
2431e5c89e4eSSatish Balay .  opt - option name
2432a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
2433a2b725a8SWilliam Gropp -  ntext - number of choices
2434e5c89e4eSSatish Balay 
2435d8d19677SJose E. Roman    Output Parameters:
24362efd9cb1SBarry Smith +  value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2437811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2438e5c89e4eSSatish Balay 
2439e5c89e4eSSatish Balay    Level: intermediate
2440e5c89e4eSSatish Balay 
244195452b02SPatrick Sanan    Notes:
244295452b02SPatrick Sanan     If the user does not supply the option value is NOT changed. Thus
24432efd9cb1SBarry Smith      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
24442efd9cb1SBarry Smith 
2445811af0c4SBarry Smith    See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList`
2446e5c89e4eSSatish Balay 
2447db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2448db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2449db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2450c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2451db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2452db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2453e5c89e4eSSatish Balay @*/
2454d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEList(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscInt ntext, PetscInt *value, PetscBool *set)
2455d71ae5a4SJacob Faibussowitsch {
245658b0ac4eSStefano Zampini   size_t    alen, len = 0, tlen = 0;
2457e5c89e4eSSatish Balay   char     *svalue;
2458ace3abfcSBarry Smith   PetscBool aset, flg = PETSC_FALSE;
2459e5c89e4eSSatish Balay   PetscInt  i;
2460e5c89e4eSSatish Balay 
2461e5c89e4eSSatish Balay   PetscFunctionBegin;
24622d747510SLisandro Dalcin   PetscValidCharPointer(opt, 3);
2463e5c89e4eSSatish Balay   for (i = 0; i < ntext; i++) {
24649566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(list[i], &alen));
2465e5c89e4eSSatish Balay     if (alen > len) len = alen;
246658b0ac4eSStefano Zampini     tlen += len + 1;
2467e5c89e4eSSatish Balay   }
2468e5c89e4eSSatish Balay   len += 5; /* a little extra space for user mistypes */
24699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &svalue));
24709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(options, pre, opt, svalue, len, &aset));
2471e5c89e4eSSatish Balay   if (aset) {
24729566063dSJacob Faibussowitsch     PetscCall(PetscEListFind(ntext, list, svalue, value, &flg));
247358b0ac4eSStefano Zampini     if (!flg) {
2474c6a7a370SJeremy L Thompson       char *avail;
247558b0ac4eSStefano Zampini 
24769566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(tlen, &avail));
2477c6a7a370SJeremy L Thompson       avail[0] = '\0';
247858b0ac4eSStefano Zampini       for (i = 0; i < ntext; i++) {
2479c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(avail, list[i], tlen));
2480c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(avail, " ", tlen));
248158b0ac4eSStefano Zampini       }
24829566063dSJacob Faibussowitsch       PetscCall(PetscStrtolower(avail));
248398921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown option %s for -%s%s. Available options: %s", svalue, pre ? pre : "", opt + 1, avail);
248458b0ac4eSStefano Zampini     }
2485fbedd5e0SJed Brown     if (set) *set = PETSC_TRUE;
2486a297a907SKarl Rupp   } else if (set) *set = PETSC_FALSE;
24879566063dSJacob Faibussowitsch   PetscCall(PetscFree(svalue));
24883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2489e5c89e4eSSatish Balay }
2490e5c89e4eSSatish Balay 
2491e5c89e4eSSatish Balay /*@C
2492e5c89e4eSSatish Balay    PetscOptionsGetEnum - Gets the enum value for a particular option in the database.
2493e5c89e4eSSatish Balay 
2494e5c89e4eSSatish Balay    Not Collective
2495e5c89e4eSSatish Balay 
2496e5c89e4eSSatish Balay    Input Parameters:
249720f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
249820f4b53cSBarry Smith .  pre - option prefix or `NULL`
2499e5c89e4eSSatish Balay .  opt - option name
25006b867d5aSJose E. Roman -  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2501e5c89e4eSSatish Balay 
2502d8d19677SJose E. Roman    Output Parameters:
2503e5c89e4eSSatish Balay +  value - the  value to return
2504811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2505e5c89e4eSSatish Balay 
2506e5c89e4eSSatish Balay    Level: beginner
2507e5c89e4eSSatish Balay 
250895452b02SPatrick Sanan    Notes:
250995452b02SPatrick Sanan     If the user does not supply the option value is NOT changed. Thus
25102efd9cb1SBarry Smith      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2511e5c89e4eSSatish Balay 
2512811af0c4SBarry Smith           List is usually something like `PCASMTypes` or some other predefined list of enum names
2513e5c89e4eSSatish Balay 
2514db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2515db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2516db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
2517db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2518c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2519db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2520db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2521e5c89e4eSSatish Balay @*/
2522d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEnum(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscEnum *value, PetscBool *set)
2523d71ae5a4SJacob Faibussowitsch {
252469a24498SJed Brown   PetscInt  ntext = 0, tval;
2525ace3abfcSBarry Smith   PetscBool fset;
2526e5c89e4eSSatish Balay 
2527e5c89e4eSSatish Balay   PetscFunctionBegin;
25282d747510SLisandro Dalcin   PetscValidCharPointer(opt, 3);
2529ad540459SPierre Jolivet   while (list[ntext++]) PetscCheck(ntext <= 50, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument appears to be wrong or have more than 50 entries");
253008401ef6SPierre Jolivet   PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix");
2531e5c89e4eSSatish Balay   ntext -= 3;
25329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEList(options, pre, opt, list, ntext, &tval, &fset));
253369a24498SJed Brown   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2534809ceb46SBarry Smith   if (fset) *value = (PetscEnum)tval;
2535809ceb46SBarry Smith   if (set) *set = fset;
25363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2537e5c89e4eSSatish Balay }
2538e5c89e4eSSatish Balay 
2539e5c89e4eSSatish Balay /*@C
25402d747510SLisandro Dalcin    PetscOptionsGetInt - Gets the integer value for a particular option in the database.
2541e5c89e4eSSatish Balay 
2542e5c89e4eSSatish Balay    Not Collective
2543e5c89e4eSSatish Balay 
2544e5c89e4eSSatish Balay    Input Parameters:
254520f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
254620f4b53cSBarry Smith .  pre - the string to prepend to the name or `NULL`
2547e5c89e4eSSatish Balay -  name - the option one is seeking
2548e5c89e4eSSatish Balay 
2549d8d19677SJose E. Roman    Output Parameters:
25502d747510SLisandro Dalcin +  ivalue - the integer value to return
2551811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2552e5c89e4eSSatish Balay 
2553e5c89e4eSSatish Balay    Level: beginner
2554e5c89e4eSSatish Balay 
2555e5c89e4eSSatish Balay    Notes:
25562d747510SLisandro Dalcin    If the user does not supply the option ivalue is NOT changed. Thus
25572efd9cb1SBarry Smith    you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
25585c07ccb8SBarry Smith 
2559db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2560db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2561db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
2562db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2563c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2564db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2565db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2566e5c89e4eSSatish Balay @*/
2567d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetInt(PetscOptions options, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
2568d71ae5a4SJacob Faibussowitsch {
25692d747510SLisandro Dalcin   const char *value;
25702d747510SLisandro Dalcin   PetscBool   flag;
2571e5c89e4eSSatish Balay 
2572e5c89e4eSSatish Balay   PetscFunctionBegin;
25732d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
25742d747510SLisandro Dalcin   PetscValidIntPointer(ivalue, 4);
25759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2576e5c89e4eSSatish Balay   if (flag) {
257734a9cc2cSBarry Smith     if (!value) {
25782d747510SLisandro Dalcin       if (set) *set = PETSC_FALSE;
257934a9cc2cSBarry Smith     } else {
25802d747510SLisandro Dalcin       if (set) *set = PETSC_TRUE;
25819566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToInt(value, ivalue));
2582e5c89e4eSSatish Balay     }
2583e5c89e4eSSatish Balay   } else {
258496ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2585e5c89e4eSSatish Balay   }
25863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2587e5c89e4eSSatish Balay }
2588e5c89e4eSSatish Balay 
2589e2446a98SMatthew Knepley /*@C
2590e5c89e4eSSatish Balay    PetscOptionsGetReal - Gets the double precision value for a particular
2591e5c89e4eSSatish Balay    option in the database.
2592e5c89e4eSSatish Balay 
2593e5c89e4eSSatish Balay    Not Collective
2594e5c89e4eSSatish Balay 
2595e5c89e4eSSatish Balay    Input Parameters:
259620f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
259720f4b53cSBarry Smith .  pre - string to prepend to each name or `NULL`
2598e5c89e4eSSatish Balay -  name - the option one is seeking
2599e5c89e4eSSatish Balay 
2600d8d19677SJose E. Roman    Output Parameters:
2601e5c89e4eSSatish Balay +  dvalue - the double value to return
2602811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, `PETSC_FALSE` if not found
2603e5c89e4eSSatish Balay 
260420f4b53cSBarry Smith    Level: beginner
260520f4b53cSBarry Smith 
2606811af0c4SBarry Smith    Note:
260795452b02SPatrick Sanan     If the user does not supply the option dvalue is NOT changed. Thus
26082efd9cb1SBarry Smith      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2609e4974155SBarry Smith 
2610db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2611c2e3fba1SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2612db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2613c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2614db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2615db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2616e5c89e4eSSatish Balay @*/
2617d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetReal(PetscOptions options, const char pre[], const char name[], PetscReal *dvalue, PetscBool *set)
2618d71ae5a4SJacob Faibussowitsch {
26192d747510SLisandro Dalcin   const char *value;
2620ace3abfcSBarry Smith   PetscBool   flag;
2621e5c89e4eSSatish Balay 
2622e5c89e4eSSatish Balay   PetscFunctionBegin;
26232d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
26242d747510SLisandro Dalcin   PetscValidRealPointer(dvalue, 4);
26259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2626e5c89e4eSSatish Balay   if (flag) {
2627a297a907SKarl Rupp     if (!value) {
2628a297a907SKarl Rupp       if (set) *set = PETSC_FALSE;
2629a297a907SKarl Rupp     } else {
2630a297a907SKarl Rupp       if (set) *set = PETSC_TRUE;
26319566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToReal(value, dvalue));
2632a297a907SKarl Rupp     }
2633e5c89e4eSSatish Balay   } else {
263496ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2635e5c89e4eSSatish Balay   }
26363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2637e5c89e4eSSatish Balay }
2638e5c89e4eSSatish Balay 
2639e5c89e4eSSatish Balay /*@C
2640e5c89e4eSSatish Balay    PetscOptionsGetScalar - Gets the scalar value for a particular
2641e5c89e4eSSatish Balay    option in the database.
2642e5c89e4eSSatish Balay 
2643e5c89e4eSSatish Balay    Not Collective
2644e5c89e4eSSatish Balay 
2645e5c89e4eSSatish Balay    Input Parameters:
264620f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
264720f4b53cSBarry Smith .  pre - string to prepend to each name or `NULL`
2648e5c89e4eSSatish Balay -  name - the option one is seeking
2649e5c89e4eSSatish Balay 
2650d8d19677SJose E. Roman    Output Parameters:
2651e5c89e4eSSatish Balay +  dvalue - the double value to return
2652811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2653e5c89e4eSSatish Balay 
2654e5c89e4eSSatish Balay    Level: beginner
2655e5c89e4eSSatish Balay 
2656e5c89e4eSSatish Balay    Usage:
2657eb4ae41dSBarry Smith    A complex number 2+3i must be specified with NO spaces
2658e5c89e4eSSatish Balay 
2659811af0c4SBarry Smith    Note:
266095452b02SPatrick Sanan     If the user does not supply the option dvalue is NOT changed. Thus
26612efd9cb1SBarry Smith      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2662e4974155SBarry Smith 
2663db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2664db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2665db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2666c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2667db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2668db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2669e5c89e4eSSatish Balay @*/
2670d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetScalar(PetscOptions options, const char pre[], const char name[], PetscScalar *dvalue, PetscBool *set)
2671d71ae5a4SJacob Faibussowitsch {
26722d747510SLisandro Dalcin   const char *value;
2673ace3abfcSBarry Smith   PetscBool   flag;
2674e5c89e4eSSatish Balay 
2675e5c89e4eSSatish Balay   PetscFunctionBegin;
26762d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
26772d747510SLisandro Dalcin   PetscValidScalarPointer(dvalue, 4);
26789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2679e5c89e4eSSatish Balay   if (flag) {
2680e5c89e4eSSatish Balay     if (!value) {
268196ef3cdfSSatish Balay       if (set) *set = PETSC_FALSE;
2682e5c89e4eSSatish Balay     } else {
2683e5c89e4eSSatish Balay #if !defined(PETSC_USE_COMPLEX)
26849566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToReal(value, dvalue));
2685e5c89e4eSSatish Balay #else
26869566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToScalar(value, dvalue));
2687e5c89e4eSSatish Balay #endif
268896ef3cdfSSatish Balay       if (set) *set = PETSC_TRUE;
2689e5c89e4eSSatish Balay     }
2690e5c89e4eSSatish Balay   } else { /* flag */
269196ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2692e5c89e4eSSatish Balay   }
26933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2694e5c89e4eSSatish Balay }
2695e5c89e4eSSatish Balay 
2696e5c89e4eSSatish Balay /*@C
2697e5c89e4eSSatish Balay    PetscOptionsGetString - Gets the string value for a particular option in
2698e5c89e4eSSatish Balay    the database.
2699e5c89e4eSSatish Balay 
2700e5c89e4eSSatish Balay    Not Collective
2701e5c89e4eSSatish Balay 
2702e5c89e4eSSatish Balay    Input Parameters:
270320f4b53cSBarry Smith +  options - options database, use `NULL` for default global database
270420f4b53cSBarry Smith .  pre - string to prepend to name or `NULL`
2705e5c89e4eSSatish Balay .  name - the option one is seeking
2706bcbf2dc5SJed Brown -  len - maximum length of the string including null termination
2707e5c89e4eSSatish Balay 
2708e5c89e4eSSatish Balay    Output Parameters:
2709e5c89e4eSSatish Balay +  string - location to copy string
2710811af0c4SBarry Smith -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2711e5c89e4eSSatish Balay 
2712e5c89e4eSSatish Balay    Level: beginner
2713e5c89e4eSSatish Balay 
271420f4b53cSBarry Smith    Note:
271520f4b53cSBarry Smith     if the option is given but no string is provided then an empty string is returned and set is given the value of `PETSC_TRUE`
271620f4b53cSBarry Smith 
271720f4b53cSBarry Smith            If the user does not use the option then the string is not changed. Thus
271820f4b53cSBarry Smith            you should ALWAYS initialize the string if you access it without first checking if the set flag is true.
271920f4b53cSBarry Smith 
272020f4b53cSBarry Smith       Even if the user provided no string (for example -optionname -someotheroption) the flag is set to PETSC_TRUE (and the string is fulled with nulls).
272120f4b53cSBarry Smith 
2722e5c89e4eSSatish Balay    Fortran Note:
2723e5c89e4eSSatish Balay    The Fortran interface is slightly different from the C/C++
2724e5c89e4eSSatish Balay    interface (len is not used).  Sample usage in Fortran follows
2725e5c89e4eSSatish Balay .vb
2726e5c89e4eSSatish Balay       character *20    string
272793e6ba5cSBarry Smith       PetscErrorCode   ierr
272893e6ba5cSBarry Smith       PetscBool        set
27291b266c99SBarry Smith       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2730e5c89e4eSSatish Balay .ve
2731e5c89e4eSSatish Balay 
2732db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2733db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2734db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2735c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2736db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2737db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2738e5c89e4eSSatish Balay @*/
2739d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetString(PetscOptions options, const char pre[], const char name[], char string[], size_t len, PetscBool *set)
2740d71ae5a4SJacob Faibussowitsch {
27412d747510SLisandro Dalcin   const char *value;
2742ace3abfcSBarry Smith   PetscBool   flag;
2743e5c89e4eSSatish Balay 
2744e5c89e4eSSatish Balay   PetscFunctionBegin;
27452d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
27462d747510SLisandro Dalcin   PetscValidCharPointer(string, 4);
27479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2748e5c89e4eSSatish Balay   if (!flag) {
274996ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2750e5c89e4eSSatish Balay   } else {
275196ef3cdfSSatish Balay     if (set) *set = PETSC_TRUE;
27529566063dSJacob Faibussowitsch     if (value) PetscCall(PetscStrncpy(string, value, len));
27539566063dSJacob Faibussowitsch     else PetscCall(PetscArrayzero(string, len));
2754e5c89e4eSSatish Balay   }
27553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2756e5c89e4eSSatish Balay }
2757e5c89e4eSSatish Balay 
2758d71ae5a4SJacob Faibussowitsch char *PetscOptionsGetStringMatlab(PetscOptions options, const char pre[], const char name[])
2759d71ae5a4SJacob Faibussowitsch {
27602d747510SLisandro Dalcin   const char *value;
276114ce751eSBarry Smith   PetscBool   flag;
276214ce751eSBarry Smith 
276314ce751eSBarry Smith   PetscFunctionBegin;
276439a651e2SJacob Faibussowitsch   if (PetscOptionsFindPair(options, pre, name, &value, &flag)) PetscFunctionReturn(NULL);
27652d747510SLisandro Dalcin   if (flag) PetscFunctionReturn((char *)value);
276639a651e2SJacob Faibussowitsch   PetscFunctionReturn(NULL);
276714ce751eSBarry Smith }
276814ce751eSBarry Smith 
27692d747510SLisandro Dalcin /*@C
27702d747510SLisandro Dalcin   PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2771f1a722f8SMatthew G. Knepley   option in the database.  The values must be separated with commas with no intervening spaces.
27722d747510SLisandro Dalcin 
27732d747510SLisandro Dalcin   Not Collective
27742d747510SLisandro Dalcin 
27752d747510SLisandro Dalcin   Input Parameters:
277620f4b53cSBarry Smith + options - options database, use `NULL` for default global database
277720f4b53cSBarry Smith . pre - string to prepend to each name or `NULL`
27786b867d5aSJose E. Roman - name - the option one is seeking
27796b867d5aSJose E. Roman 
2780d8d19677SJose E. Roman   Output Parameters:
27812d747510SLisandro Dalcin + dvalue - the integer values to return
2782f1a722f8SMatthew G. Knepley . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2783811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
27842d747510SLisandro Dalcin 
27852d747510SLisandro Dalcin   Level: beginner
27862d747510SLisandro Dalcin 
2787811af0c4SBarry Smith   Note:
278820f4b53cSBarry Smith   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`. FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
27892d747510SLisandro Dalcin 
2790db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2791db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2792db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2793c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2794db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2795db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
27962d747510SLisandro Dalcin @*/
2797d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options, const char pre[], const char name[], PetscBool dvalue[], PetscInt *nmax, PetscBool *set)
2798d71ae5a4SJacob Faibussowitsch {
27992d747510SLisandro Dalcin   const char *svalue;
28002d747510SLisandro Dalcin   char       *value;
28012d747510SLisandro Dalcin   PetscInt    n = 0;
28022d747510SLisandro Dalcin   PetscBool   flag;
28032d747510SLisandro Dalcin   PetscToken  token;
28042d747510SLisandro Dalcin 
28052d747510SLisandro Dalcin   PetscFunctionBegin;
28062d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
2807064a246eSJacob Faibussowitsch   PetscValidBoolPointer(dvalue, 4);
28082d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 5);
28092d747510SLisandro Dalcin 
28109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
28119371c9d4SSatish Balay   if (!flag || !svalue) {
28129371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
28139371c9d4SSatish Balay     *nmax = 0;
28143ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
28159371c9d4SSatish Balay   }
28162d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
28179566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
28189566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
28192d747510SLisandro Dalcin   while (value && n < *nmax) {
28209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(value, dvalue));
28219566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
28222d747510SLisandro Dalcin     dvalue++;
28232d747510SLisandro Dalcin     n++;
28242d747510SLisandro Dalcin   }
28259566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
28262d747510SLisandro Dalcin   *nmax = n;
28273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28282d747510SLisandro Dalcin }
28292d747510SLisandro Dalcin 
28302d747510SLisandro Dalcin /*@C
28312d747510SLisandro Dalcin   PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.
28322d747510SLisandro Dalcin 
28332d747510SLisandro Dalcin   Not Collective
28342d747510SLisandro Dalcin 
28352d747510SLisandro Dalcin   Input Parameters:
283620f4b53cSBarry Smith + options - options database, use `NULL` for default global database
283720f4b53cSBarry Smith . pre - option prefix or `NULL`
28382d747510SLisandro Dalcin . name - option name
28396b867d5aSJose E. Roman - list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
28406b867d5aSJose E. Roman 
28412d747510SLisandro Dalcin   Output Parameters:
28422d747510SLisandro Dalcin + ivalue - the  enum values to return
2843f1a722f8SMatthew G. Knepley . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2844811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
28452d747510SLisandro Dalcin 
28462d747510SLisandro Dalcin   Level: beginner
28472d747510SLisandro Dalcin 
28482d747510SLisandro Dalcin   Notes:
28492d747510SLisandro Dalcin   The array must be passed as a comma separated list.
28502d747510SLisandro Dalcin 
28512d747510SLisandro Dalcin   There must be no intervening spaces between the values.
28522d747510SLisandro Dalcin 
2853811af0c4SBarry Smith   list is usually something like `PCASMTypes` or some other predefined list of enum names.
28542d747510SLisandro Dalcin 
2855db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2856db781477SPatrick Sanan           `PetscOptionsGetEnum()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2857db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, `PetscOptionsName()`,
2858c2e3fba1SPatrick Sanan           `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, `PetscOptionsStringArray()`, `PetscOptionsRealArray()`,
2859db781477SPatrick Sanan           `PetscOptionsScalar()`, `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2860db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
28612d747510SLisandro Dalcin @*/
2862d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options, const char pre[], const char name[], const char *const *list, PetscEnum ivalue[], PetscInt *nmax, PetscBool *set)
2863d71ae5a4SJacob Faibussowitsch {
28642d747510SLisandro Dalcin   const char *svalue;
28652d747510SLisandro Dalcin   char       *value;
28662d747510SLisandro Dalcin   PetscInt    n = 0;
28672d747510SLisandro Dalcin   PetscEnum   evalue;
28682d747510SLisandro Dalcin   PetscBool   flag;
28692d747510SLisandro Dalcin   PetscToken  token;
28702d747510SLisandro Dalcin 
28712d747510SLisandro Dalcin   PetscFunctionBegin;
28722d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
28732d747510SLisandro Dalcin   PetscValidPointer(list, 4);
28742d747510SLisandro Dalcin   PetscValidPointer(ivalue, 5);
28752d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 6);
28762d747510SLisandro Dalcin 
28779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
28789371c9d4SSatish Balay   if (!flag || !svalue) {
28799371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
28809371c9d4SSatish Balay     *nmax = 0;
28813ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
28829371c9d4SSatish Balay   }
28832d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
28849566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
28859566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
28862d747510SLisandro Dalcin   while (value && n < *nmax) {
28879566063dSJacob Faibussowitsch     PetscCall(PetscEnumFind(list, value, &evalue, &flag));
288828b400f6SJacob Faibussowitsch     PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown enum value '%s' for -%s%s", svalue, pre ? pre : "", name + 1);
28892d747510SLisandro Dalcin     ivalue[n++] = evalue;
28909566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
28912d747510SLisandro Dalcin   }
28929566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
28932d747510SLisandro Dalcin   *nmax = n;
28943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28952d747510SLisandro Dalcin }
28962d747510SLisandro Dalcin 
28972d747510SLisandro Dalcin /*@C
2898f1a722f8SMatthew G. Knepley   PetscOptionsGetIntArray - Gets an array of integer values for a particular option in the database.
28992d747510SLisandro Dalcin 
29002d747510SLisandro Dalcin   Not Collective
29012d747510SLisandro Dalcin 
29022d747510SLisandro Dalcin   Input Parameters:
290320f4b53cSBarry Smith + options - options database, use `NULL` for default global database
290420f4b53cSBarry Smith . pre - string to prepend to each name or `NULL`
29056b867d5aSJose E. Roman - name - the option one is seeking
29066b867d5aSJose E. Roman 
2907d8d19677SJose E. Roman   Output Parameters:
29082d747510SLisandro Dalcin + ivalue - the integer values to return
2909f1a722f8SMatthew G. Knepley . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2910811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
29112d747510SLisandro Dalcin 
29122d747510SLisandro Dalcin   Level: beginner
29132d747510SLisandro Dalcin 
29142d747510SLisandro Dalcin   Notes:
29152d747510SLisandro Dalcin   The array can be passed as
2916811af0c4SBarry Smith +  a comma separated list -                                 0,1,2,3,4,5,6,7
2917811af0c4SBarry Smith .  a range (start\-end+1) -                                 0-8
2918811af0c4SBarry Smith .  a range with given increment (start\-end+1:inc) -        0-7:2
2919811af0c4SBarry Smith -  a combination of values and ranges separated by commas - 0,1-8,8-15:2
29202d747510SLisandro Dalcin 
29212d747510SLisandro Dalcin   There must be no intervening spaces between the values.
29222d747510SLisandro Dalcin 
2923db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2924db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2925db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2926c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2927db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2928db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
29292d747510SLisandro Dalcin @*/
2930d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetIntArray(PetscOptions options, const char pre[], const char name[], PetscInt ivalue[], PetscInt *nmax, PetscBool *set)
2931d71ae5a4SJacob Faibussowitsch {
29322d747510SLisandro Dalcin   const char *svalue;
29332d747510SLisandro Dalcin   char       *value;
29342d747510SLisandro Dalcin   PetscInt    n = 0, i, j, start, end, inc, nvalues;
29352d747510SLisandro Dalcin   size_t      len;
29362d747510SLisandro Dalcin   PetscBool   flag, foundrange;
29372d747510SLisandro Dalcin   PetscToken  token;
29382d747510SLisandro Dalcin 
29392d747510SLisandro Dalcin   PetscFunctionBegin;
29402d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
29412d747510SLisandro Dalcin   PetscValidIntPointer(ivalue, 4);
29422d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 5);
29432d747510SLisandro Dalcin 
29449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
29459371c9d4SSatish Balay   if (!flag || !svalue) {
29469371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
29479371c9d4SSatish Balay     *nmax = 0;
29483ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
29499371c9d4SSatish Balay   }
29502d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
29519566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
29529566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
29532d747510SLisandro Dalcin   while (value && n < *nmax) {
29542d747510SLisandro Dalcin     /* look for form  d-D where d and D are integers */
29552d747510SLisandro Dalcin     foundrange = PETSC_FALSE;
29569566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(value, &len));
29572d747510SLisandro Dalcin     if (value[0] == '-') i = 2;
29582d747510SLisandro Dalcin     else i = 1;
29592d747510SLisandro Dalcin     for (; i < (int)len; i++) {
29602d747510SLisandro Dalcin       if (value[i] == '-') {
2961cc73adaaSBarry Smith         PetscCheck(i != (int)len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value);
29622d747510SLisandro Dalcin         value[i] = 0;
29632d747510SLisandro Dalcin 
29649566063dSJacob Faibussowitsch         PetscCall(PetscOptionsStringToInt(value, &start));
29652d747510SLisandro Dalcin         inc = 1;
29662d747510SLisandro Dalcin         j   = i + 1;
29672d747510SLisandro Dalcin         for (; j < (int)len; j++) {
29682d747510SLisandro Dalcin           if (value[j] == ':') {
29692d747510SLisandro Dalcin             value[j] = 0;
29702d747510SLisandro Dalcin 
29719566063dSJacob Faibussowitsch             PetscCall(PetscOptionsStringToInt(value + j + 1, &inc));
297208401ef6SPierre Jolivet             PetscCheck(inc > 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry,%s cannot have negative increment", n, value + j + 1);
29732d747510SLisandro Dalcin             break;
29742d747510SLisandro Dalcin           }
29752d747510SLisandro Dalcin         }
29769566063dSJacob Faibussowitsch         PetscCall(PetscOptionsStringToInt(value + i + 1, &end));
297708401ef6SPierre Jolivet         PetscCheck(end > start, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry, %s-%s cannot have decreasing list", n, value, value + i + 1);
29782d747510SLisandro Dalcin         nvalues = (end - start) / inc + (end - start) % inc;
2979cc73adaaSBarry Smith         PetscCheck(n + nvalues <= *nmax, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry, not enough space left in array (%" PetscInt_FMT ") to contain entire range from %" PetscInt_FMT " to %" PetscInt_FMT, n, *nmax - n, start, end);
29802d747510SLisandro Dalcin         for (; start < end; start += inc) {
29819371c9d4SSatish Balay           *ivalue = start;
29829371c9d4SSatish Balay           ivalue++;
29839371c9d4SSatish Balay           n++;
29842d747510SLisandro Dalcin         }
29852d747510SLisandro Dalcin         foundrange = PETSC_TRUE;
29862d747510SLisandro Dalcin         break;
29872d747510SLisandro Dalcin       }
29882d747510SLisandro Dalcin     }
29892d747510SLisandro Dalcin     if (!foundrange) {
29909566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToInt(value, ivalue));
29912d747510SLisandro Dalcin       ivalue++;
29922d747510SLisandro Dalcin       n++;
29932d747510SLisandro Dalcin     }
29949566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
29952d747510SLisandro Dalcin   }
29969566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
29972d747510SLisandro Dalcin   *nmax = n;
29983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29992d747510SLisandro Dalcin }
30002d747510SLisandro Dalcin 
30012d747510SLisandro Dalcin /*@C
30022d747510SLisandro Dalcin   PetscOptionsGetRealArray - Gets an array of double precision values for a
3003f1a722f8SMatthew G. Knepley   particular option in the database.  The values must be separated with commas with no intervening spaces.
30042d747510SLisandro Dalcin 
30052d747510SLisandro Dalcin   Not Collective
30062d747510SLisandro Dalcin 
30072d747510SLisandro Dalcin   Input Parameters:
300820f4b53cSBarry Smith + options - options database, use `NULL` for default global database
300920f4b53cSBarry Smith . pre - string to prepend to each name or `NULL`
30106b867d5aSJose E. Roman - name - the option one is seeking
30116b867d5aSJose E. Roman 
30122d747510SLisandro Dalcin   Output Parameters:
30132d747510SLisandro Dalcin + dvalue - the double values to return
3014f1a722f8SMatthew G. Knepley . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
3015811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
30162d747510SLisandro Dalcin 
30172d747510SLisandro Dalcin   Level: beginner
30182d747510SLisandro Dalcin 
3019db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3020db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3021db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3022c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3023db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3024db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
30252d747510SLisandro Dalcin @*/
3026d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetRealArray(PetscOptions options, const char pre[], const char name[], PetscReal dvalue[], PetscInt *nmax, PetscBool *set)
3027d71ae5a4SJacob Faibussowitsch {
30282d747510SLisandro Dalcin   const char *svalue;
30292d747510SLisandro Dalcin   char       *value;
30302d747510SLisandro Dalcin   PetscInt    n = 0;
30312d747510SLisandro Dalcin   PetscBool   flag;
30322d747510SLisandro Dalcin   PetscToken  token;
30332d747510SLisandro Dalcin 
30342d747510SLisandro Dalcin   PetscFunctionBegin;
30352d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
30362d747510SLisandro Dalcin   PetscValidRealPointer(dvalue, 4);
30372d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 5);
30382d747510SLisandro Dalcin 
30399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
30409371c9d4SSatish Balay   if (!flag || !svalue) {
30419371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
30429371c9d4SSatish Balay     *nmax = 0;
30433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
30449371c9d4SSatish Balay   }
30452d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
30469566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
30479566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
30482d747510SLisandro Dalcin   while (value && n < *nmax) {
30499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToReal(value, dvalue++));
30509566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
30512d747510SLisandro Dalcin     n++;
30522d747510SLisandro Dalcin   }
30539566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
30542d747510SLisandro Dalcin   *nmax = n;
30553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30562d747510SLisandro Dalcin }
30572d747510SLisandro Dalcin 
30582d747510SLisandro Dalcin /*@C
30592d747510SLisandro Dalcin   PetscOptionsGetScalarArray - Gets an array of scalars for a
3060f1a722f8SMatthew G. Knepley   particular option in the database.  The values must be separated with commas with no intervening spaces.
30612d747510SLisandro Dalcin 
30622d747510SLisandro Dalcin   Not Collective
30632d747510SLisandro Dalcin 
30642d747510SLisandro Dalcin   Input Parameters:
306520f4b53cSBarry Smith + options - options database, use `NULL` for default global database
306620f4b53cSBarry Smith . pre - string to prepend to each name or `NULL`
30676b867d5aSJose E. Roman - name - the option one is seeking
30686b867d5aSJose E. Roman 
30692d747510SLisandro Dalcin   Output Parameters:
30702d747510SLisandro Dalcin + dvalue - the scalar values to return
3071f1a722f8SMatthew G. Knepley . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
3072811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
30732d747510SLisandro Dalcin 
30742d747510SLisandro Dalcin   Level: beginner
30752d747510SLisandro Dalcin 
3076db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3077db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3078db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3079c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3080db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3081db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
30822d747510SLisandro Dalcin @*/
3083d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options, const char pre[], const char name[], PetscScalar dvalue[], PetscInt *nmax, PetscBool *set)
3084d71ae5a4SJacob Faibussowitsch {
30852d747510SLisandro Dalcin   const char *svalue;
30862d747510SLisandro Dalcin   char       *value;
30872d747510SLisandro Dalcin   PetscInt    n = 0;
30882d747510SLisandro Dalcin   PetscBool   flag;
30892d747510SLisandro Dalcin   PetscToken  token;
30902d747510SLisandro Dalcin 
30912d747510SLisandro Dalcin   PetscFunctionBegin;
30922d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
3093064a246eSJacob Faibussowitsch   PetscValidScalarPointer(dvalue, 4);
30942d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 5);
30952d747510SLisandro Dalcin 
30969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
30979371c9d4SSatish Balay   if (!flag || !svalue) {
30989371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
30999371c9d4SSatish Balay     *nmax = 0;
31003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
31019371c9d4SSatish Balay   }
31022d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
31039566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
31049566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
31052d747510SLisandro Dalcin   while (value && n < *nmax) {
31069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToScalar(value, dvalue++));
31079566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
31082d747510SLisandro Dalcin     n++;
31092d747510SLisandro Dalcin   }
31109566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
31112d747510SLisandro Dalcin   *nmax = n;
31123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31132d747510SLisandro Dalcin }
311414ce751eSBarry Smith 
3115e5c89e4eSSatish Balay /*@C
3116e5c89e4eSSatish Balay   PetscOptionsGetStringArray - Gets an array of string values for a particular
3117f1a722f8SMatthew G. Knepley   option in the database. The values must be separated with commas with no intervening spaces.
3118e5c89e4eSSatish Balay 
3119cf53795eSBarry Smith   Not Collective; No Fortran Support
3120e5c89e4eSSatish Balay 
3121e5c89e4eSSatish Balay   Input Parameters:
312220f4b53cSBarry Smith + options - options database, use `NULL` for default global database
312320f4b53cSBarry Smith . pre - string to prepend to name or `NULL`
31246b867d5aSJose E. Roman - name - the option one is seeking
31256b867d5aSJose E. Roman 
3126e7b76fa7SPatrick Sanan   Output Parameters:
3127e5c89e4eSSatish Balay + strings - location to copy strings
3128f1a722f8SMatthew G. Knepley . nmax - On input maximum number of strings, on output the actual number of strings found
3129811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
3130e5c89e4eSSatish Balay 
3131e5c89e4eSSatish Balay   Level: beginner
3132e5c89e4eSSatish Balay 
3133e5c89e4eSSatish Balay   Notes:
3134e7b76fa7SPatrick Sanan   The nmax parameter is used for both input and output.
3135e7b76fa7SPatrick Sanan 
3136e5c89e4eSSatish Balay   The user should pass in an array of pointers to char, to hold all the
3137e5c89e4eSSatish Balay   strings returned by this function.
3138e5c89e4eSSatish Balay 
3139e5c89e4eSSatish Balay   The user is responsible for deallocating the strings that are
3140cf53795eSBarry Smith   returned.
3141e5c89e4eSSatish Balay 
3142db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
3143db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3144db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3145c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3146db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3147db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
3148e5c89e4eSSatish Balay @*/
3149d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetStringArray(PetscOptions options, const char pre[], const char name[], char *strings[], PetscInt *nmax, PetscBool *set)
3150d71ae5a4SJacob Faibussowitsch {
31512d747510SLisandro Dalcin   const char *svalue;
3152e5c89e4eSSatish Balay   char       *value;
31532d747510SLisandro Dalcin   PetscInt    n = 0;
3154ace3abfcSBarry Smith   PetscBool   flag;
31559c9d3cfdSBarry Smith   PetscToken  token;
3156e5c89e4eSSatish Balay 
3157e5c89e4eSSatish Balay   PetscFunctionBegin;
31582d747510SLisandro Dalcin   PetscValidCharPointer(name, 3);
31592d747510SLisandro Dalcin   PetscValidPointer(strings, 4);
31602d747510SLisandro Dalcin   PetscValidIntPointer(nmax, 5);
3161e5c89e4eSSatish Balay 
31629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
31639371c9d4SSatish Balay   if (!flag || !svalue) {
31649371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
31659371c9d4SSatish Balay     *nmax = 0;
31663ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
31679371c9d4SSatish Balay   }
31682d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
31699566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
31709566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
31712d747510SLisandro Dalcin   while (value && n < *nmax) {
31729566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(value, &strings[n]));
31739566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
3174e5c89e4eSSatish Balay     n++;
3175e5c89e4eSSatish Balay   }
31769566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
3177e5c89e4eSSatish Balay   *nmax = n;
31783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3179e5c89e4eSSatish Balay }
318006824ed3SPatrick Sanan 
318106824ed3SPatrick Sanan /*@C
31824ead3382SBarry Smith    PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with `newname`
318306824ed3SPatrick Sanan 
318406824ed3SPatrick Sanan    Prints a deprecation warning, unless an option is supplied to suppress.
318506824ed3SPatrick Sanan 
31861c9f3c13SBarry Smith    Logically Collective
318706824ed3SPatrick Sanan 
318806824ed3SPatrick Sanan    Input Parameters:
318920f4b53cSBarry Smith +  pre - string to prepend to name or `NULL`
319006824ed3SPatrick Sanan .  oldname - the old, deprecated option
319120f4b53cSBarry Smith .  newname - the new option, or `NULL` if option is purely removed
31929f3a6782SPatrick Sanan .  version - a string describing the version of first deprecation, e.g. "3.9"
319320f4b53cSBarry Smith -  info - additional information string, or `NULL`.
319406824ed3SPatrick Sanan 
3195811af0c4SBarry Smith    Options Database Key:
319606824ed3SPatrick Sanan . -options_suppress_deprecated_warnings - do not print deprecation warnings
319706824ed3SPatrick Sanan 
319820f4b53cSBarry Smith    Level: developer
319920f4b53cSBarry Smith 
320006824ed3SPatrick Sanan    Notes:
32014ead3382SBarry Smith    If `newname` is provided then the options database will automatically check the database for `oldname`.
32024ead3382SBarry Smith 
32034ead3382SBarry Smith    The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
32044ead3382SBarry Smith    new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
32054ead3382SBarry Smith    See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
32064ead3382SBarry Smith 
3207811af0c4SBarry Smith    Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
320835cb6cd3SPierre Jolivet    Only the process of rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or
3209811af0c4SBarry Smith    `PetscObjectOptionsBegin()` prints the information
3210b40114eaSPatrick Sanan    If newname is provided, the old option is replaced. Otherwise, it remains
3211b40114eaSPatrick Sanan    in the options database.
32129f3a6782SPatrick Sanan    If an option is not replaced, the info argument should be used to advise the user
32139f3a6782SPatrick Sanan    on how to proceed.
32149f3a6782SPatrick Sanan    There is a limit on the length of the warning printed, so very long strings
32159f3a6782SPatrick Sanan    provided as info may be truncated.
321606824ed3SPatrick Sanan 
3217db781477SPatrick Sanan .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
321806824ed3SPatrick Sanan @*/
3219d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject, const char oldname[], const char newname[], const char version[], const char info[])
3220d71ae5a4SJacob Faibussowitsch {
322106824ed3SPatrick Sanan   PetscBool         found, quiet;
322206824ed3SPatrick Sanan   const char       *value;
322306824ed3SPatrick Sanan   const char *const quietopt = "-options_suppress_deprecated_warnings";
32249f3a6782SPatrick Sanan   char              msg[4096];
3225b0bdc838SStefano Zampini   char             *prefix  = NULL;
3226b0bdc838SStefano Zampini   PetscOptions      options = NULL;
3227b0bdc838SStefano Zampini   MPI_Comm          comm    = PETSC_COMM_SELF;
322806824ed3SPatrick Sanan 
322906824ed3SPatrick Sanan   PetscFunctionBegin;
323006824ed3SPatrick Sanan   PetscValidCharPointer(oldname, 2);
323106824ed3SPatrick Sanan   PetscValidCharPointer(version, 4);
3232b0bdc838SStefano Zampini   if (PetscOptionsObject) {
3233b0bdc838SStefano Zampini     prefix  = PetscOptionsObject->prefix;
3234b0bdc838SStefano Zampini     options = PetscOptionsObject->options;
3235b0bdc838SStefano Zampini     comm    = PetscOptionsObject->comm;
3236b0bdc838SStefano Zampini   }
32379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, prefix, oldname, &value, &found));
323806824ed3SPatrick Sanan   if (found) {
323906824ed3SPatrick Sanan     if (newname) {
32401baa6e33SBarry Smith       if (prefix) PetscCall(PetscOptionsPrefixPush(options, prefix));
32419566063dSJacob Faibussowitsch       PetscCall(PetscOptionsSetValue(options, newname, value));
32421baa6e33SBarry Smith       if (prefix) PetscCall(PetscOptionsPrefixPop(options));
32439566063dSJacob Faibussowitsch       PetscCall(PetscOptionsClearValue(options, oldname));
3244b40114eaSPatrick Sanan     }
324506824ed3SPatrick Sanan     quiet = PETSC_FALSE;
32469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(options, NULL, quietopt, &quiet, NULL));
324706824ed3SPatrick Sanan     if (!quiet) {
3248c6a7a370SJeremy L Thompson       PetscCall(PetscStrncpy(msg, "** PETSc DEPRECATION WARNING ** : the option ", sizeof(msg)));
3249c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, oldname, sizeof(msg)));
3250c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, " is deprecated as of version ", sizeof(msg)));
3251c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, version, sizeof(msg)));
32524bd3d7f8SBarry Smith       PetscCall(PetscStrlcat(msg, " and will be removed in a future release.\n", sizeof(msg)));
325306824ed3SPatrick Sanan       if (newname) {
32544bd3d7f8SBarry Smith         PetscCall(PetscStrlcat(msg, "   Use the option ", sizeof(msg)));
3255c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, newname, sizeof(msg)));
3256c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, " instead.", sizeof(msg)));
325706824ed3SPatrick Sanan       }
32589f3a6782SPatrick Sanan       if (info) {
3259c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, " ", sizeof(msg)));
3260c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, info, sizeof(msg)));
32619f3a6782SPatrick Sanan       }
3262c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, " (Silence this warning with ", sizeof(msg)));
3263c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, quietopt, sizeof(msg)));
3264c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, ")\n", sizeof(msg)));
32659566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%s", msg));
326606824ed3SPatrick Sanan     }
326706824ed3SPatrick Sanan   }
32683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326906824ed3SPatrick Sanan }
3270