xref: /petsc/src/sys/objects/options.c (revision 432b765a2304277672d65f4442022404e0cbdb25)
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;
1279355ec05SMatthew G. Knepley   if (options->monitorFromOptions) PetscCall(PetscOptionsMonitorDefault(name, value, source, NULL));
1289355ec05SMatthew G. Knepley   for (PetscInt i = 0; i < options->numbermonitors; i++) PetscCall((*options->monitor[i])(name, value, source, options->monitorcontext[i]));
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130e5c89e4eSSatish Balay }
131e5c89e4eSSatish Balay 
1322d747510SLisandro Dalcin /*@
1332d747510SLisandro Dalcin   PetscOptionsCreate - Creates an empty options database.
134e5c89e4eSSatish Balay 
13520f4b53cSBarry Smith   Logically Collective
1361c9f3c13SBarry Smith 
137e5c89e4eSSatish Balay   Output Parameter:
1382d747510SLisandro Dalcin . options - Options database object
139e5c89e4eSSatish Balay 
140e5c89e4eSSatish Balay   Level: advanced
141e5c89e4eSSatish Balay 
142811af0c4SBarry Smith   Note:
143811af0c4SBarry Smith   Though PETSc has a concept of multiple options database the current code uses a single default `PetscOptions` object
144811af0c4SBarry Smith 
145811af0c4SBarry Smith   Developer Notes:
146811af0c4SBarry Smith   We may want eventually to pass a `MPI_Comm` to determine the ownership of the object
147811af0c4SBarry Smith 
148811af0c4SBarry Smith   This object never got developed after being introduced, it is not clear that supporting multiple `PetscOptions` objects is useful
1491c9f3c13SBarry Smith 
150db781477SPatrick Sanan .seealso: `PetscOptionsDestroy()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
151e5c89e4eSSatish Balay @*/
152d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsCreate(PetscOptions *options)
153d71ae5a4SJacob Faibussowitsch {
15439a651e2SJacob Faibussowitsch   PetscFunctionBegin;
1554f572ea9SToby Isaac   PetscAssertPointer(options, 1);
1562d747510SLisandro Dalcin   *options = (PetscOptions)calloc(1, sizeof(**options));
15739a651e2SJacob Faibussowitsch   PetscCheck(*options, PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate the options database");
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1592d747510SLisandro Dalcin }
1602d747510SLisandro Dalcin 
1612d747510SLisandro Dalcin /*@
1622d747510SLisandro Dalcin   PetscOptionsDestroy - Destroys an option database.
1632d747510SLisandro Dalcin 
16420f4b53cSBarry Smith   Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
1651c9f3c13SBarry Smith 
1662d747510SLisandro Dalcin   Input Parameter:
167811af0c4SBarry Smith . options - the `PetscOptions` object
1682d747510SLisandro Dalcin 
1693de2bfdfSBarry Smith   Level: advanced
1702d747510SLisandro Dalcin 
171aec76313SJacob Faibussowitsch .seealso: `PetscOptionsInsert()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsSetValue()`
1722d747510SLisandro Dalcin @*/
173d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
174d71ae5a4SJacob Faibussowitsch {
175362febeeSStefano Zampini   PetscFunctionBegin;
1764f572ea9SToby Isaac   PetscAssertPointer(options, 1);
1773ba16761SJacob Faibussowitsch   if (!*options) PetscFunctionReturn(PETSC_SUCCESS);
1785f80ce2aSJacob 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()");
1799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsClear(*options));
1802d747510SLisandro Dalcin   /* XXX what about monitors ? */
1812800570dSLisandro Dalcin   free(*options);
1822d747510SLisandro Dalcin   *options = NULL;
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184e5c89e4eSSatish Balay }
185e5c89e4eSSatish Balay 
1862d747510SLisandro Dalcin /*
1872d747510SLisandro Dalcin     PetscOptionsCreateDefault - Creates the default global options database
1882d747510SLisandro Dalcin */
189d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsCreateDefault(void)
190d71ae5a4SJacob Faibussowitsch {
19139a651e2SJacob Faibussowitsch   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   if (PetscUnlikely(!defaultoptions)) PetscCall(PetscOptionsCreate(&defaultoptions));
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1942d747510SLisandro Dalcin }
1952d747510SLisandro Dalcin 
196b4205f0bSBarry Smith /*@
197811af0c4SBarry Smith   PetscOptionsPush - Push a new `PetscOptions` object as the default provider of options
1981c9f3c13SBarry Smith   Allows using different parts of a code to use different options databases
199b4205f0bSBarry Smith 
200b4205f0bSBarry Smith   Logically Collective
201b4205f0bSBarry Smith 
202b4205f0bSBarry Smith   Input Parameter:
203811af0c4SBarry Smith . opt - the options obtained with `PetscOptionsCreate()`
204b4205f0bSBarry Smith 
20520f4b53cSBarry Smith   Level: advanced
20620f4b53cSBarry Smith 
207b4205f0bSBarry Smith   Notes:
208811af0c4SBarry Smith   Use `PetscOptionsPop()` to return to the previous default options database
2091c9f3c13SBarry Smith 
210811af0c4SBarry Smith   The collectivity of this routine is complex; only the MPI ranks that call this routine will
2111c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
2121c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
2131c9f3c13SBarry Smith   on different ranks.
214b4205f0bSBarry Smith 
215aec76313SJacob Faibussowitsch   Developer Notes:
216811af0c4SBarry Smith   Though this functionality has been provided it has never been used in PETSc and might be removed.
217811af0c4SBarry Smith 
218db781477SPatrick Sanan .seealso: `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
219b4205f0bSBarry Smith @*/
220d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPush(PetscOptions opt)
221d71ae5a4SJacob Faibussowitsch {
222b4205f0bSBarry Smith   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsCreateDefault());
224b4205f0bSBarry Smith   opt->previous  = defaultoptions;
225b4205f0bSBarry Smith   defaultoptions = opt;
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
227b4205f0bSBarry Smith }
228b4205f0bSBarry Smith 
229b4205f0bSBarry Smith /*@
230811af0c4SBarry Smith   PetscOptionsPop - Pop the most recent `PetscOptionsPush()` to return to the previous default options
231b4205f0bSBarry Smith 
23220f4b53cSBarry Smith   Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
233b4205f0bSBarry Smith 
2343de2bfdfSBarry Smith   Level: advanced
2353de2bfdfSBarry Smith 
23642747ad1SJacob Faibussowitsch .seealso: `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
237b4205f0bSBarry Smith @*/
238d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPop(void)
239d71ae5a4SJacob Faibussowitsch {
2403de2bfdfSBarry Smith   PetscOptions current = defaultoptions;
2413de2bfdfSBarry Smith 
242b4205f0bSBarry Smith   PetscFunctionBegin;
24328b400f6SJacob Faibussowitsch   PetscCheck(defaultoptions, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing default options");
24428b400f6SJacob Faibussowitsch   PetscCheck(defaultoptions->previous, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscOptionsPop() called too many times");
245b4205f0bSBarry Smith   defaultoptions    = defaultoptions->previous;
2463de2bfdfSBarry Smith   current->previous = NULL;
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
248b4205f0bSBarry Smith }
249b4205f0bSBarry Smith 
2502d747510SLisandro Dalcin /*
2512d747510SLisandro Dalcin     PetscOptionsDestroyDefault - Destroys the default global options database
2522d747510SLisandro Dalcin */
253d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDestroyDefault(void)
254d71ae5a4SJacob Faibussowitsch {
25539a651e2SJacob Faibussowitsch   PetscFunctionBegin;
2563ba16761SJacob Faibussowitsch   if (!defaultoptions) PetscFunctionReturn(PETSC_SUCCESS);
2573de2bfdfSBarry Smith   /* Destroy any options that the user forgot to pop */
2583de2bfdfSBarry Smith   while (defaultoptions->previous) {
25939a651e2SJacob Faibussowitsch     PetscOptions tmp = defaultoptions;
26039a651e2SJacob Faibussowitsch 
2619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsPop());
2629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsDestroy(&tmp));
2633de2bfdfSBarry Smith   }
2649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDestroy(&defaultoptions));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266e5c89e4eSSatish Balay }
267e5c89e4eSSatish Balay 
26894ef8ddeSSatish Balay /*@C
2697cd08cecSJed Brown   PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.
2703fc1eb6aSBarry Smith 
27120f4b53cSBarry Smith   Not Collective
2721c9f3c13SBarry Smith 
2733fc1eb6aSBarry Smith   Input Parameter:
2742d747510SLisandro Dalcin . key - string to check if valid
2753fc1eb6aSBarry Smith 
2763fc1eb6aSBarry Smith   Output Parameter:
277811af0c4SBarry Smith . valid - `PETSC_TRUE` if a valid key
2783fc1eb6aSBarry Smith 
279f6680f47SSatish Balay   Level: intermediate
28010450e9eSJacob Faibussowitsch 
28110450e9eSJacob Faibussowitsch .seealso: `PetscOptionsCreate()`, `PetscOptionsInsert()`
2823fc1eb6aSBarry Smith @*/
283d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsValidKey(const char key[], PetscBool *valid)
284d71ae5a4SJacob Faibussowitsch {
285f603b5e9SToby Isaac   char               *ptr;
28627304958SStefano Zampini   PETSC_UNUSED double d;
2877c5db45bSBarry Smith 
28896fc60bcSBarry Smith   PetscFunctionBegin;
2894f572ea9SToby Isaac   if (key) PetscAssertPointer(key, 1);
2904f572ea9SToby Isaac   PetscAssertPointer(valid, 2);
2912d747510SLisandro Dalcin   *valid = PETSC_FALSE;
2923ba16761SJacob Faibussowitsch   if (!key) PetscFunctionReturn(PETSC_SUCCESS);
2933ba16761SJacob Faibussowitsch   if (key[0] != '-') PetscFunctionReturn(PETSC_SUCCESS);
2942d747510SLisandro Dalcin   if (key[1] == '-') key++;
2953ba16761SJacob Faibussowitsch   if (!isalpha((int)key[1])) PetscFunctionReturn(PETSC_SUCCESS);
29627304958SStefano Zampini   d = strtod(key, &ptr);
2973ba16761SJacob Faibussowitsch   if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) PetscFunctionReturn(PETSC_SUCCESS);
2982d747510SLisandro Dalcin   *valid = PETSC_TRUE;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30096fc60bcSBarry Smith }
30196fc60bcSBarry Smith 
30210c654e6SJacob Faibussowitsch static PetscErrorCode PetscOptionsInsertString_Private(PetscOptions options, const char in_str[], PetscOptionSource source)
303d71ae5a4SJacob Faibussowitsch {
304d06005cbSLisandro Dalcin   char      *first, *second;
3059c9d3cfdSBarry Smith   PetscToken token;
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(in_str, ' ', &token));
3099566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &first));
31096fc60bcSBarry Smith   while (first) {
311d06005cbSLisandro Dalcin     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
31210c654e6SJacob Faibussowitsch 
3139566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_file", &isfile));
3149566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_file_yaml", &isfileyaml));
3159566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-options_string_yaml", &isstringyaml));
3169566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-prefix_push", &ispush));
3179566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(first, "-prefix_pop", &ispop));
3189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(first, &key));
319d06005cbSLisandro Dalcin     if (!key) {
3209566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
321d06005cbSLisandro Dalcin     } else if (isfile) {
3229566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
32310c654e6SJacob Faibussowitsch       PetscCall(PetscOptionsInsertFile(PETSC_COMM_SELF, options, second, PETSC_TRUE));
3249566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
325d06005cbSLisandro Dalcin     } else if (isfileyaml) {
3269566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
32710c654e6SJacob Faibussowitsch       PetscCall(PetscOptionsInsertFileYAML(PETSC_COMM_SELF, options, second, PETSC_TRUE));
3289566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
329d06005cbSLisandro Dalcin     } else if (isstringyaml) {
3309566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3319355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, second, source));
3329566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
333d06005cbSLisandro Dalcin     } else if (ispush) {
3349566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3359566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPush(options, second));
3369566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
3379db968c8SJed Brown     } else if (ispop) {
3389566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPop(options));
3399566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &first));
340d06005cbSLisandro Dalcin     } else {
3419566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &second));
3429566063dSJacob Faibussowitsch       PetscCall(PetscOptionsValidKey(second, &key));
34396fc60bcSBarry Smith       if (!key) {
3449355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, first, second, NULL, source));
3459566063dSJacob Faibussowitsch         PetscCall(PetscTokenFind(token, &first));
34696fc60bcSBarry Smith       } else {
3479355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, first, NULL, NULL, source));
34896fc60bcSBarry Smith         first = second;
34996fc60bcSBarry Smith       }
350e5c89e4eSSatish Balay     }
351e5c89e4eSSatish Balay   }
3529566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354e5c89e4eSSatish Balay }
355e5c89e4eSSatish Balay 
3569355ec05SMatthew G. Knepley /*@C
3579355ec05SMatthew G. Knepley   PetscOptionsInsertString - Inserts options into the database from a string
3589355ec05SMatthew G. Knepley 
3599355ec05SMatthew G. Knepley   Logically Collective
3609355ec05SMatthew G. Knepley 
3619355ec05SMatthew G. Knepley   Input Parameters:
3629355ec05SMatthew G. Knepley + options - options object
3639355ec05SMatthew G. Knepley - in_str  - string that contains options separated by blanks
3649355ec05SMatthew G. Knepley 
3659355ec05SMatthew G. Knepley   Level: intermediate
3669355ec05SMatthew G. Knepley 
36720f4b53cSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
3689355ec05SMatthew G. Knepley   have the affect of these options. If some processes that create objects call this routine and others do
3699355ec05SMatthew G. Knepley   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
3709355ec05SMatthew G. Knepley   on different ranks.
3719355ec05SMatthew G. Knepley 
3729355ec05SMatthew G. Knepley    Contributed by Boyana Norris
3739355ec05SMatthew G. Knepley 
3749355ec05SMatthew G. Knepley .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
3759355ec05SMatthew G. Knepley           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3769355ec05SMatthew G. Knepley           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3779355ec05SMatthew G. Knepley           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3789355ec05SMatthew G. Knepley           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3799355ec05SMatthew G. Knepley           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsInsertFile()`
3809355ec05SMatthew G. Knepley @*/
3819355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsInsertString(PetscOptions options, const char in_str[])
3829355ec05SMatthew G. Knepley {
3839355ec05SMatthew G. Knepley   PetscFunctionBegin;
3849355ec05SMatthew G. Knepley   PetscCall(PetscOptionsInsertString_Private(options, in_str, PETSC_OPT_CODE));
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3869355ec05SMatthew G. Knepley }
3879355ec05SMatthew G. Knepley 
3883fc1eb6aSBarry Smith /*
3893fc1eb6aSBarry Smith     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
3903fc1eb6aSBarry Smith */
391d71ae5a4SJacob Faibussowitsch static char *Petscgetline(FILE *f)
392d71ae5a4SJacob Faibussowitsch {
3935fa91da5SBarry Smith   size_t size = 0;
3945fa91da5SBarry Smith   size_t len  = 0;
3955fa91da5SBarry Smith   size_t last = 0;
3960298fd71SBarry Smith   char  *buf  = NULL;
3975fa91da5SBarry Smith 
39802c9f0b5SLisandro Dalcin   if (feof(f)) return NULL;
3995fa91da5SBarry Smith   do {
4005fa91da5SBarry Smith     size += 1024;                             /* BUFSIZ is defined as "the optimal read size for this platform" */
4016e0c8459SSatish Balay     buf = (char *)realloc((void *)buf, size); /* realloc(NULL,n) is the same as malloc(n) */
4025fa91da5SBarry Smith     /* Actually do the read. Note that fgets puts a terminal '\0' on the
4035fa91da5SBarry Smith     end of the string, so we make sure we overwrite this */
404e86f3e45SDave May     if (!fgets(buf + len, 1024, f)) buf[len] = 0;
4053ba16761SJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, PetscStrlen(buf, &len));
4065fa91da5SBarry Smith     last = len - 1;
4075fa91da5SBarry Smith   } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
40808ac41f7SSatish Balay   if (len) return buf;
4095fa91da5SBarry Smith   free(buf);
41002c9f0b5SLisandro Dalcin   return NULL;
4115fa91da5SBarry Smith }
4125fa91da5SBarry Smith 
413d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsFilename(MPI_Comm comm, const char file[], char filename[PETSC_MAX_PATH_LEN], PetscBool *yaml)
414d71ae5a4SJacob Faibussowitsch {
415be10d61cSLisandro Dalcin   char fname[PETSC_MAX_PATH_LEN + 8], path[PETSC_MAX_PATH_LEN + 8], *tail;
416e5c89e4eSSatish Balay 
417be10d61cSLisandro Dalcin   PetscFunctionBegin;
418362febeeSStefano Zampini   *yaml = PETSC_FALSE;
4199566063dSJacob Faibussowitsch   PetscCall(PetscStrreplace(comm, file, fname, sizeof(fname)));
4209566063dSJacob Faibussowitsch   PetscCall(PetscFixFilename(fname, path));
4219566063dSJacob Faibussowitsch   PetscCall(PetscStrendswith(path, ":yaml", yaml));
422be10d61cSLisandro Dalcin   if (*yaml) {
4239566063dSJacob Faibussowitsch     PetscCall(PetscStrrchr(path, ':', &tail));
424be10d61cSLisandro Dalcin     tail[-1] = 0; /* remove ":yaml" suffix from path */
425be10d61cSLisandro Dalcin   }
4269566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename, path, PETSC_MAX_PATH_LEN));
427a1d2f846SLisandro Dalcin   /* check for standard YAML and JSON filename extensions */
4289566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yaml", yaml));
4299566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yml", yaml));
4309566063dSJacob Faibussowitsch   if (!*yaml) PetscCall(PetscStrendswith(filename, ".json", yaml));
431a1d2f846SLisandro Dalcin   if (!*yaml) { /* check file contents */
432a1d2f846SLisandro Dalcin     PetscMPIInt rank;
4339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
434dd400576SPatrick Sanan     if (rank == 0) {
435a1d2f846SLisandro Dalcin       FILE *fh = fopen(filename, "r");
436a1d2f846SLisandro Dalcin       if (fh) {
437a1d2f846SLisandro Dalcin         char buf[6] = "";
438a1d2f846SLisandro Dalcin         if (fread(buf, 1, 6, fh) > 0) {
4399566063dSJacob Faibussowitsch           PetscCall(PetscStrncmp(buf, "%YAML ", 6, yaml));          /* check for '%YAML' tag */
4409566063dSJacob Faibussowitsch           if (!*yaml) PetscCall(PetscStrncmp(buf, "---", 3, yaml)); /* check for document start */
441a1d2f846SLisandro Dalcin         }
442a1d2f846SLisandro Dalcin         (void)fclose(fh);
443a1d2f846SLisandro Dalcin       }
444a1d2f846SLisandro Dalcin     }
4459566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(yaml, 1, MPIU_BOOL, 0, comm));
446a1d2f846SLisandro Dalcin   }
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
448be10d61cSLisandro Dalcin }
449e5c89e4eSSatish Balay 
450d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsInsertFilePetsc(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
451d71ae5a4SJacob Faibussowitsch {
4528c0b561eSLisandro Dalcin   char       *string, *vstring = NULL, *astring = NULL, *packed = NULL;
4537fb43599SVaclav Hapla   char       *tokens[4];
45413e3f751SJed Brown   size_t      i, len, bytes;
455e5c89e4eSSatish Balay   FILE       *fd;
4567fb43599SVaclav Hapla   PetscToken  token = NULL;
457ed9cf6e9SBarry Smith   int         err;
458bbcf679cSJacob Faibussowitsch   char       *cmatch = NULL;
459581bbe83SVaclav Hapla   const char  cmt    = '#';
4609210b8eaSVaclav Hapla   PetscInt    line   = 1;
4613a018368SJed Brown   PetscMPIInt rank, cnt = 0, acnt = 0, counts[2];
4629210b8eaSVaclav Hapla   PetscBool   isdir, alias = PETSC_FALSE, valid;
463e5c89e4eSSatish Balay 
464e5c89e4eSSatish Balay   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(tokens, sizeof(tokens)));
4669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
467dd400576SPatrick Sanan   if (rank == 0) {
4688c0b561eSLisandro Dalcin     char fpath[PETSC_MAX_PATH_LEN];
4698c0b561eSLisandro Dalcin     char fname[PETSC_MAX_PATH_LEN];
47005c7dedfSBarry Smith 
4719566063dSJacob Faibussowitsch     PetscCall(PetscStrreplace(PETSC_COMM_SELF, file, fpath, sizeof(fpath)));
4729566063dSJacob Faibussowitsch     PetscCall(PetscFixFilename(fpath, fname));
4738c0b561eSLisandro Dalcin 
474e5c89e4eSSatish Balay     fd = fopen(fname, "r");
4759566063dSJacob Faibussowitsch     PetscCall(PetscTestDirectory(fname, 'r', &isdir));
47608401ef6SPierre Jolivet     PetscCheck(!isdir || !require, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified options file %s is a directory", fname);
477ad38b122SPatrick Sanan     if (fd && !isdir) {
4783a018368SJed Brown       PetscSegBuffer vseg, aseg;
4799566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferCreate(1, 4000, &vseg));
4809566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferCreate(1, 2000, &aseg));
4813a018368SJed Brown 
4829b754dc9SBarry Smith       /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
4839566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Opened options file %s\n", file));
484e24ecc5dSJed Brown 
4855fa91da5SBarry Smith       while ((string = Petscgetline(fd))) {
4864704e885SBarry Smith         /* eliminate comments from each line */
4879566063dSJacob Faibussowitsch         PetscCall(PetscStrchr(string, cmt, &cmatch));
48890f79514SSatish Balay         if (cmatch) *cmatch = 0;
4899566063dSJacob Faibussowitsch         PetscCall(PetscStrlen(string, &len));
4905981331cSSatish Balay         /* replace tabs, ^M, \n with " " */
491e5c89e4eSSatish Balay         for (i = 0; i < len; i++) {
492ad540459SPierre Jolivet           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') string[i] = ' ';
493e5c89e4eSSatish Balay         }
4949566063dSJacob Faibussowitsch         PetscCall(PetscTokenCreate(string, ' ', &token));
4959566063dSJacob Faibussowitsch         PetscCall(PetscTokenFind(token, &tokens[0]));
4967fb43599SVaclav Hapla         if (!tokens[0]) {
49702b0d46eSSatish Balay           goto destroy;
4987fb43599SVaclav Hapla         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
4999566063dSJacob Faibussowitsch           PetscCall(PetscTokenFind(token, &tokens[0]));
50090f79514SSatish Balay         }
50148a46eb9SPierre Jolivet         for (i = 1; i < 4; i++) PetscCall(PetscTokenFind(token, &tokens[i]));
5027fb43599SVaclav Hapla         if (!tokens[0]) {
5032662f744SSatish Balay           goto destroy;
5047fb43599SVaclav Hapla         } else if (tokens[0][0] == '-') {
5059566063dSJacob Faibussowitsch           PetscCall(PetscOptionsValidKey(tokens[0], &valid));
50628b400f6SJacob 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]);
5079566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(tokens[0], &len));
5089566063dSJacob Faibussowitsch           PetscCall(PetscSegBufferGet(vseg, len + 1, &vstring));
5099566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(vstring, tokens[0], len));
510e24ecc5dSJed Brown           vstring[len] = ' ';
5117fb43599SVaclav Hapla           if (tokens[1]) {
5129566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
51328b400f6SJacob 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]);
5149566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[1], &len));
5159566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(vseg, len + 3, &vstring));
516e24ecc5dSJed Brown             vstring[0] = '"';
5179566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(vstring + 1, tokens[1], len));
518e24ecc5dSJed Brown             vstring[len + 1] = '"';
519e24ecc5dSJed Brown             vstring[len + 2] = ' ';
52009192fe3SBarry Smith           }
52190f79514SSatish Balay         } else {
5229566063dSJacob Faibussowitsch           PetscCall(PetscStrcasecmp(tokens[0], "alias", &alias));
5239210b8eaSVaclav Hapla           if (alias) {
5249566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
52528b400f6SJacob 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]);
52608401ef6SPierre 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]);
5279566063dSJacob Faibussowitsch             PetscCall(PetscOptionsValidKey(tokens[2], &valid));
52828b400f6SJacob 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]);
5299566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[1], &len));
5309566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
5319566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(astring, tokens[1], len));
532e24ecc5dSJed Brown             astring[len] = ' ';
533e24ecc5dSJed Brown 
5349566063dSJacob Faibussowitsch             PetscCall(PetscStrlen(tokens[2], &len));
5359566063dSJacob Faibussowitsch             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
5369566063dSJacob Faibussowitsch             PetscCall(PetscArraycpy(astring, tokens[2], len));
537e24ecc5dSJed Brown             astring[len] = ' ';
53898921bdaSJacob 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]);
5399210b8eaSVaclav Hapla         }
5409210b8eaSVaclav Hapla         {
5419210b8eaSVaclav Hapla           const char *extraToken = alias ? tokens[3] : tokens[2];
54228b400f6SJacob Faibussowitsch           PetscCheck(!extraToken, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": extra token %s", fname, line, extraToken);
543e5c89e4eSSatish Balay         }
54402b0d46eSSatish Balay       destroy:
5454b40f50bSBarry Smith         free(string);
5469566063dSJacob Faibussowitsch         PetscCall(PetscTokenDestroy(&token));
5479210b8eaSVaclav Hapla         alias = PETSC_FALSE;
5489210b8eaSVaclav Hapla         line++;
549e5c89e4eSSatish Balay       }
550ed9cf6e9SBarry Smith       err = fclose(fd);
55128b400f6SJacob Faibussowitsch       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file %s", fname);
5529566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGetSize(aseg, &bytes)); /* size without null termination */
5539566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(bytes, &acnt));
5549566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(aseg, 1, &astring));
555e24ecc5dSJed Brown       astring[0] = 0;
5569566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGetSize(vseg, &bytes)); /* size without null termination */
5579566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(bytes, &cnt));
5589566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(vseg, 1, &vstring));
559e24ecc5dSJed Brown       vstring[0] = 0;
5609566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
5619566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferExtractTo(aseg, packed));
5629566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferExtractTo(vseg, packed + acnt + 1));
5639566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferDestroy(&aseg));
5649566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferDestroy(&vseg));
56528b400f6SJacob Faibussowitsch     } else PetscCheck(!require, PETSC_COMM_SELF, PETSC_ERR_USER, "Unable to open options file %s", fname);
5669b754dc9SBarry Smith   }
56705c7dedfSBarry Smith 
5683a018368SJed Brown   counts[0] = acnt;
5693a018368SJed Brown   counts[1] = cnt;
5704201f521SBarry Smith   err       = MPI_Bcast(counts, 2, MPI_INT, 0, comm);
57128b400f6SJacob 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/");
5723a018368SJed Brown   acnt = counts[0];
5733a018368SJed Brown   cnt  = counts[1];
57448a46eb9SPierre Jolivet   if (rank) PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
5753a018368SJed Brown   if (acnt || cnt) {
5769566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(packed, 2 + acnt + cnt, MPI_CHAR, 0, comm));
5773a018368SJed Brown     astring = packed;
5783a018368SJed Brown     vstring = packed + acnt + 1;
5793a018368SJed Brown   }
5803a018368SJed Brown 
5819b754dc9SBarry Smith   if (acnt) {
5829566063dSJacob Faibussowitsch     PetscCall(PetscTokenCreate(astring, ' ', &token));
5839566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &tokens[0]));
5847fb43599SVaclav Hapla     while (tokens[0]) {
5859566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &tokens[1]));
5869566063dSJacob Faibussowitsch       PetscCall(PetscOptionsSetAlias(options, tokens[0], tokens[1]));
5879566063dSJacob Faibussowitsch       PetscCall(PetscTokenFind(token, &tokens[0]));
5889b754dc9SBarry Smith     }
5899566063dSJacob Faibussowitsch     PetscCall(PetscTokenDestroy(&token));
5909b754dc9SBarry Smith   }
5919b754dc9SBarry Smith 
5929355ec05SMatthew G. Knepley   if (cnt) PetscCall(PetscOptionsInsertString_Private(options, vstring, PETSC_OPT_FILE));
5939566063dSJacob Faibussowitsch   PetscCall(PetscFree(packed));
5943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
595e5c89e4eSSatish Balay }
596e5c89e4eSSatish Balay 
597d06005cbSLisandro Dalcin /*@C
598be10d61cSLisandro Dalcin   PetscOptionsInsertFile - Inserts options into the database from a file.
599be10d61cSLisandro Dalcin 
600be10d61cSLisandro Dalcin   Collective
601be10d61cSLisandro Dalcin 
602d8d19677SJose E. Roman   Input Parameters:
603811af0c4SBarry Smith + comm    - the processes that will share the options (usually `PETSC_COMM_WORLD`)
60420f4b53cSBarry Smith . options - options database, use `NULL` for default global database
605be10d61cSLisandro Dalcin . file    - name of file,
606be10d61cSLisandro Dalcin            ".yml" and ".yaml" filename extensions are inserted as YAML options,
607be10d61cSLisandro Dalcin            append ":yaml" to filename to force YAML options.
608811af0c4SBarry Smith - require - if `PETSC_TRUE` will generate an error if the file does not exist
609be10d61cSLisandro Dalcin 
61020f4b53cSBarry Smith   Level: developer
61120f4b53cSBarry Smith 
612be10d61cSLisandro Dalcin   Notes:
613be10d61cSLisandro Dalcin   Use  # for lines that are comments and which should be ignored.
614811af0c4SBarry Smith   Usually, instead of using this command, one should list the file name in the call to `PetscInitialize()`, this insures that certain options
61521532e8aSBarry 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
61621532e8aSBarry Smith   calls to `XXXSetFromOptions()`, it should not be used for options listed under PetscInitialize().
61721532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes in comm will
61821532e8aSBarry Smith   have the effect of these options. If some processes that create objects call this routine and others do
619be10d61cSLisandro Dalcin   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
620be10d61cSLisandro Dalcin   on different ranks.
621be10d61cSLisandro Dalcin 
622db781477SPatrick Sanan .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
623db781477SPatrick Sanan           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
624db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
625c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
626db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
627db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
628be10d61cSLisandro Dalcin @*/
629d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
630d71ae5a4SJacob Faibussowitsch {
631be10d61cSLisandro Dalcin   char      filename[PETSC_MAX_PATH_LEN];
632be10d61cSLisandro Dalcin   PetscBool yaml;
633be10d61cSLisandro Dalcin 
634be10d61cSLisandro Dalcin   PetscFunctionBegin;
6359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFilename(comm, file, filename, &yaml));
636be10d61cSLisandro Dalcin   if (yaml) {
6379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFileYAML(comm, options, filename, require));
638be10d61cSLisandro Dalcin   } else {
6399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFilePetsc(comm, options, filename, require));
640be10d61cSLisandro Dalcin   }
6413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
642be10d61cSLisandro Dalcin }
643be10d61cSLisandro Dalcin 
644be10d61cSLisandro Dalcin /*@C
645d06005cbSLisandro Dalcin   PetscOptionsInsertArgs - Inserts options into the database from a array of strings
646d06005cbSLisandro Dalcin 
647d06005cbSLisandro Dalcin   Logically Collective
648d06005cbSLisandro Dalcin 
649d8d19677SJose E. Roman   Input Parameters:
650d06005cbSLisandro Dalcin + options - options object
6516aad120cSJose E. Roman . argc    - the array length
652d06005cbSLisandro Dalcin - args    - the string array
653d06005cbSLisandro Dalcin 
654d06005cbSLisandro Dalcin   Level: intermediate
655d06005cbSLisandro Dalcin 
656db781477SPatrick Sanan .seealso: `PetscOptions`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`
657d06005cbSLisandro Dalcin @*/
658d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsertArgs(PetscOptions options, int argc, char *args[])
659d71ae5a4SJacob Faibussowitsch {
660d06005cbSLisandro Dalcin   MPI_Comm     comm  = PETSC_COMM_WORLD;
661d06005cbSLisandro Dalcin   int          left  = PetscMax(argc, 0);
662d06005cbSLisandro Dalcin   char *const *eargs = args;
66385079163SJed Brown 
66485079163SJed Brown   PetscFunctionBegin;
66585079163SJed Brown   while (left) {
666d06005cbSLisandro Dalcin     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
6679566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_file", &isfile));
6689566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_file_yaml", &isfileyaml));
6699566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-options_string_yaml", &isstringyaml));
6709566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_push", &ispush));
6719566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_pop", &ispop));
6729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(eargs[0], &key));
673093de6efSBarry Smith     if (!key) {
6749371c9d4SSatish Balay       eargs++;
6759371c9d4SSatish Balay       left--;
676d06005cbSLisandro Dalcin     } else if (isfile) {
677cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file filename option");
6789566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInsertFile(comm, options, eargs[1], PETSC_TRUE));
6799371c9d4SSatish Balay       eargs += 2;
6809371c9d4SSatish Balay       left -= 2;
681d06005cbSLisandro Dalcin     } else if (isfileyaml) {
682cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file_yaml filename option");
6839566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInsertFileYAML(comm, options, eargs[1], PETSC_TRUE));
6849371c9d4SSatish Balay       eargs += 2;
6859371c9d4SSatish Balay       left -= 2;
686d06005cbSLisandro Dalcin     } else if (isstringyaml) {
687cc73adaaSBarry Smith       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing string for -options_string_yaml string option");
6889355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, eargs[1], PETSC_OPT_CODE));
6899371c9d4SSatish Balay       eargs += 2;
6909371c9d4SSatish Balay       left -= 2;
691d06005cbSLisandro Dalcin     } else if (ispush) {
69208401ef6SPierre Jolivet       PetscCheck(left > 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option");
693cc73adaaSBarry Smith       PetscCheck(eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option (prefixes cannot start with '-')");
6949566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPush(options, eargs[1]));
6959371c9d4SSatish Balay       eargs += 2;
6969371c9d4SSatish Balay       left -= 2;
697d06005cbSLisandro Dalcin     } else if (ispop) {
6989566063dSJacob Faibussowitsch       PetscCall(PetscOptionsPrefixPop(options));
6999371c9d4SSatish Balay       eargs++;
7009371c9d4SSatish Balay       left--;
7017935c3d8SJed Brown     } else {
7027935c3d8SJed Brown       PetscBool nextiskey = PETSC_FALSE;
7039566063dSJacob Faibussowitsch       if (left >= 2) PetscCall(PetscOptionsValidKey(eargs[1], &nextiskey));
70498b6bf53SJed Brown       if (left < 2 || nextiskey) {
7059355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], NULL, NULL, PETSC_OPT_COMMAND_LINE));
7069371c9d4SSatish Balay         eargs++;
7079371c9d4SSatish Balay         left--;
70885079163SJed Brown       } else {
7099355ec05SMatthew G. Knepley         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], eargs[1], NULL, PETSC_OPT_COMMAND_LINE));
7109371c9d4SSatish Balay         eargs += 2;
7119371c9d4SSatish Balay         left -= 2;
71285079163SJed Brown       }
71385079163SJed Brown     }
7147935c3d8SJed Brown   }
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71685079163SJed Brown }
71785079163SJed Brown 
71810c654e6SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt, const char *val[], const PetscBool set[], PetscBool *flg)
719d71ae5a4SJacob Faibussowitsch {
720c5b5d8d5SVaclav Hapla   PetscFunctionBegin;
721c5b5d8d5SVaclav Hapla   if (set[opt]) {
7229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(val[opt], flg));
723c5b5d8d5SVaclav Hapla   } else *flg = PETSC_FALSE;
7243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
725c5b5d8d5SVaclav Hapla }
726c5b5d8d5SVaclav Hapla 
727660278c0SBarry Smith /* Process options with absolute precedence, these are only processed from the command line, not the environment or files */
728d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options, int argc, char *args[], PetscBool *skip_petscrc, PetscBool *skip_petscrc_set)
729d71ae5a4SJacob Faibussowitsch {
730c5b5d8d5SVaclav Hapla   const char *const *opt = precedentOptions;
731c5b5d8d5SVaclav Hapla   const size_t       n   = PO_NUM;
732c5b5d8d5SVaclav Hapla   size_t             o;
733c5b5d8d5SVaclav Hapla   int                a;
734c5b5d8d5SVaclav Hapla   const char       **val;
7350c99d500SBarry Smith   char             **cval;
736660278c0SBarry Smith   PetscBool         *set, unneeded;
737c5b5d8d5SVaclav Hapla 
738c5b5d8d5SVaclav Hapla   PetscFunctionBegin;
7390c99d500SBarry Smith   PetscCall(PetscCalloc2(n, &cval, n, &set));
7400c99d500SBarry Smith   val = (const char **)cval;
741c5b5d8d5SVaclav Hapla 
742c5b5d8d5SVaclav Hapla   /* Look for options possibly set using PetscOptionsSetValue beforehand */
74348a46eb9SPierre Jolivet   for (o = 0; o < n; o++) PetscCall(PetscOptionsFindPair(options, NULL, opt[o], &val[o], &set[o]));
744c5b5d8d5SVaclav Hapla 
745a5b23f4aSJose E. Roman   /* Loop through all args to collect last occurring value of each option */
746c5b5d8d5SVaclav Hapla   for (a = 1; a < argc; a++) {
747c5b5d8d5SVaclav Hapla     PetscBool valid, eq;
748c5b5d8d5SVaclav Hapla 
7499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(args[a], &valid));
750c5b5d8d5SVaclav Hapla     if (!valid) continue;
751c5b5d8d5SVaclav Hapla     for (o = 0; o < n; o++) {
7529566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(args[a], opt[o], &eq));
753c5b5d8d5SVaclav Hapla       if (eq) {
754c5b5d8d5SVaclav Hapla         set[o] = PETSC_TRUE;
755c5b5d8d5SVaclav Hapla         if (a == argc - 1 || !args[a + 1] || !args[a + 1][0] || args[a + 1][0] == '-') val[o] = NULL;
756c5b5d8d5SVaclav Hapla         else val[o] = args[a + 1];
757c5b5d8d5SVaclav Hapla         break;
758c5b5d8d5SVaclav Hapla       }
759c5b5d8d5SVaclav Hapla     }
760c5b5d8d5SVaclav Hapla   }
761c5b5d8d5SVaclav Hapla 
762c5b5d8d5SVaclav Hapla   /* Process flags */
7639566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro));
764d314f959SVaclav Hapla   if (options->help_intro) options->help = PETSC_TRUE;
7659566063dSJacob Faibussowitsch   else PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_HELP, val, set, &options->help));
766660278c0SBarry Smith   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_CI_ENABLE, val, set, &unneeded));
767660278c0SBarry Smith   /* need to manage PO_CI_ENABLE option before the PetscOptionsMonitor is turned on, so its setting is not monitored */
7689355ec05SMatthew G. Knepley   if (set[PO_CI_ENABLE]) PetscCall(PetscOptionsSetValue_Private(options, opt[PO_CI_ENABLE], val[PO_CI_ENABLE], &a, PETSC_OPT_COMMAND_LINE));
7699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL, val, set, &options->monitorCancel));
7709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR, val, set, &options->monitorFromOptions));
7719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC, val, set, skip_petscrc));
772c5b5d8d5SVaclav Hapla   *skip_petscrc_set = set[PO_SKIP_PETSCRC];
773c5b5d8d5SVaclav Hapla 
774c5b5d8d5SVaclav Hapla   /* Store precedent options in database and mark them as used */
775660278c0SBarry Smith   for (o = 1; o < n; o++) {
776c5b5d8d5SVaclav Hapla     if (set[o]) {
7779355ec05SMatthew G. Knepley       PetscCall(PetscOptionsSetValue_Private(options, opt[o], val[o], &a, PETSC_OPT_COMMAND_LINE));
778d06005cbSLisandro Dalcin       options->used[a] = PETSC_TRUE;
779c5b5d8d5SVaclav Hapla     }
780c5b5d8d5SVaclav Hapla   }
7810c99d500SBarry Smith   PetscCall(PetscFree2(cval, set));
782c5b5d8d5SVaclav Hapla   options->precedentProcessed = PETSC_TRUE;
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784c5b5d8d5SVaclav Hapla }
785c5b5d8d5SVaclav Hapla 
786d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options, const char name[], PetscBool *flg)
787d71ae5a4SJacob Faibussowitsch {
78839a651e2SJacob Faibussowitsch   PetscFunctionBegin;
7894f572ea9SToby Isaac   PetscAssertPointer(flg, 3);
790c5b5d8d5SVaclav Hapla   *flg = PETSC_FALSE;
791c5b5d8d5SVaclav Hapla   if (options->precedentProcessed) {
79239a651e2SJacob Faibussowitsch     for (int i = 0; i < PO_NUM; ++i) {
793c5b5d8d5SVaclav Hapla       if (!PetscOptNameCmp(precedentOptions[i], name)) {
794c5b5d8d5SVaclav Hapla         /* check if precedent option has been set already */
7959566063dSJacob Faibussowitsch         PetscCall(PetscOptionsFindPair(options, NULL, name, NULL, flg));
796c5b5d8d5SVaclav Hapla         if (*flg) break;
797c5b5d8d5SVaclav Hapla       }
798c5b5d8d5SVaclav Hapla     }
799c5b5d8d5SVaclav Hapla   }
8003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
801c5b5d8d5SVaclav Hapla }
80285079163SJed Brown 
803e5c89e4eSSatish Balay /*@C
804e5c89e4eSSatish Balay   PetscOptionsInsert - Inserts into the options database from the command line,
805e5c89e4eSSatish Balay   the environmental variable and a file.
806e5c89e4eSSatish Balay 
807811af0c4SBarry Smith   Collective on `PETSC_COMM_WORLD`
8081c9f3c13SBarry Smith 
809e5c89e4eSSatish Balay   Input Parameters:
81020f4b53cSBarry Smith + options - options database or `NULL` for the default global database
811c5929fdfSBarry Smith . argc    - count of number of command line arguments
812e5c89e4eSSatish Balay . args    - the command line arguments
813be10d61cSLisandro Dalcin - file    - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
81420f4b53cSBarry Smith           Use `NULL` or empty string to not check for code specific file.
815be10d61cSLisandro Dalcin           Also checks ~/.petscrc, .petscrc and petscrc.
816c5b5d8d5SVaclav Hapla           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
817e5c89e4eSSatish Balay 
818081c24baSBoyana Norris   Options Database Keys:
819d06005cbSLisandro Dalcin + -options_file <filename>      - read options from a file
820d06005cbSLisandro Dalcin - -options_file_yaml <filename> - read options from a YAML file
821c5b5d8d5SVaclav Hapla 
82220f4b53cSBarry Smith   Level: advanced
82320f4b53cSBarry Smith 
824811af0c4SBarry Smith   Notes:
82520f4b53cSBarry Smith   Since `PetscOptionsInsert()` is automatically called by `PetscInitialize()`,
826811af0c4SBarry Smith   the user does not typically need to call this routine. `PetscOptionsInsert()`
827811af0c4SBarry Smith   can be called several times, adding additional entries into the database.
828811af0c4SBarry Smith 
829811af0c4SBarry Smith   See `PetscInitialize()` for options related to option database monitoring.
830081c24baSBoyana Norris 
831db781477SPatrick Sanan .seealso: `PetscOptionsDestroy()`, `PetscOptionsView()`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`,
832db781477SPatrick Sanan           `PetscInitialize()`
833e5c89e4eSSatish Balay @*/
834d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInsert(PetscOptions options, int *argc, char ***args, const char file[])
835d71ae5a4SJacob Faibussowitsch {
836d06005cbSLisandro Dalcin   MPI_Comm    comm = PETSC_COMM_WORLD;
837e5c89e4eSSatish Balay   PetscMPIInt rank;
838c5b5d8d5SVaclav Hapla   PetscBool   hasArgs     = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
839c5b5d8d5SVaclav Hapla   PetscBool   skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;
840e5c89e4eSSatish Balay 
841e5c89e4eSSatish Balay   PetscFunctionBegin;
84208401ef6SPierre Jolivet   PetscCheck(!hasArgs || (args && *args), comm, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
8439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
844e5c89e4eSSatish Balay 
845c5b5d8d5SVaclav Hapla   if (!options) {
8469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsCreateDefault());
847c5b5d8d5SVaclav Hapla     options = defaultoptions;
848c5b5d8d5SVaclav Hapla   }
849c5b5d8d5SVaclav Hapla   if (hasArgs) {
850c5b5d8d5SVaclav Hapla     /* process options with absolute precedence */
8519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsProcessPrecedentFlags(options, *argc, *args, &skipPetscrc, &skipPetscrcSet));
852660278c0SBarry Smith     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci", &PetscCIEnabled, NULL));
853c5b5d8d5SVaclav Hapla   }
8544b09e917SBarry Smith   if (file && file[0]) {
8559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, file, PETSC_TRUE));
856c5b5d8d5SVaclav Hapla     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
8579566063dSJacob Faibussowitsch     if (!skipPetscrcSet) PetscCall(PetscOptionsGetBool(options, NULL, "-skip_petscrc", &skipPetscrc, NULL));
858321366bcSBarry Smith   }
859c5b5d8d5SVaclav Hapla   if (!skipPetscrc) {
860be10d61cSLisandro Dalcin     char filename[PETSC_MAX_PATH_LEN];
8619566063dSJacob Faibussowitsch     PetscCall(PetscGetHomeDirectory(filename, sizeof(filename)));
8629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(filename, (int)sizeof(filename), MPI_CHAR, 0, comm));
863c6a7a370SJeremy L Thompson     if (filename[0]) PetscCall(PetscStrlcat(filename, "/.petscrc", sizeof(filename)));
8649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, filename, PETSC_FALSE));
8659566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, ".petscrc", PETSC_FALSE));
8669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInsertFile(comm, options, "petscrc", PETSC_FALSE));
867e5c89e4eSSatish Balay   }
868e5c89e4eSSatish Balay 
8692d747510SLisandro Dalcin   /* insert environment options */
870e5c89e4eSSatish Balay   {
8712d747510SLisandro Dalcin     char  *eoptions = NULL;
872e5c89e4eSSatish Balay     size_t len      = 0;
873dd400576SPatrick Sanan     if (rank == 0) {
874e5c89e4eSSatish Balay       eoptions = (char *)getenv("PETSC_OPTIONS");
8759566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(eoptions, &len));
876e5c89e4eSSatish Balay     }
8779566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
878e5c89e4eSSatish Balay     if (len) {
8799566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
8809566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
88196fc60bcSBarry Smith       if (rank) eoptions[len] = 0;
8829355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertString_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
8839566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscFree(eoptions));
884e5c89e4eSSatish Balay     }
885e5c89e4eSSatish Balay   }
886e5c89e4eSSatish Balay 
887d06005cbSLisandro Dalcin   /* insert YAML environment options */
88856a31166SBarry Smith   {
8899fc438c3SToby Isaac     char  *eoptions = NULL;
8909fc438c3SToby Isaac     size_t len      = 0;
891dd400576SPatrick Sanan     if (rank == 0) {
8929fc438c3SToby Isaac       eoptions = (char *)getenv("PETSC_OPTIONS_YAML");
8939566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(eoptions, &len));
8949fc438c3SToby Isaac     }
8959566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
8969fc438c3SToby Isaac     if (len) {
8979566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
8989566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
8999fc438c3SToby Isaac       if (rank) eoptions[len] = 0;
9009355ec05SMatthew G. Knepley       PetscCall(PetscOptionsInsertStringYAML_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
9019566063dSJacob Faibussowitsch       if (rank) PetscCall(PetscFree(eoptions));
9029fc438c3SToby Isaac     }
9039fc438c3SToby Isaac   }
9043bcbd388SSean Farley 
905c5b5d8d5SVaclav Hapla   /* insert command line options here because they take precedence over arguments in petscrc/environment */
9069566063dSJacob Faibussowitsch   if (hasArgs) PetscCall(PetscOptionsInsertArgs(options, *argc - 1, *args + 1));
907660278c0SBarry Smith   PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci_portable_error_output", &PetscCIEnabledPortableErrorOutput, NULL));
9083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
909e5c89e4eSSatish Balay }
910e5c89e4eSSatish Balay 
911660278c0SBarry Smith /* These options are not printed with PetscOptionsView() or PetscOptionsMonitor() when PetscCIEnabled is on */
912660278c0SBarry Smith /* TODO: get the list from the test harness, do not have it hardwired here. Maybe from gmakegentest.py */
913e9a33e21SBarry Smith static const char *PetscCIOptions[] = {"malloc_debug", "malloc_dump", "malloc_test", "malloc", "nox", "nox_warning", "display", "saws_port_auto_select", "saws_port_auto_select_silent", "vecscatter_mpi1", "check_pointer_intensity", "cuda_initialize", "error_output_stdout", "use_gpu_aware_mpi", "checkfunctionlist", "fp_trap", "petsc_ci", "petsc_ci_portable_error_output", "options_left"};
914660278c0SBarry Smith 
915d71ae5a4SJacob Faibussowitsch static PetscBool PetscCIOption(const char *name)
916d71ae5a4SJacob Faibussowitsch {
917660278c0SBarry Smith   PetscInt  idx;
918660278c0SBarry Smith   PetscBool found;
919660278c0SBarry Smith 
920660278c0SBarry Smith   if (!PetscCIEnabled) return PETSC_FALSE;
9213ba16761SJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, PetscEListFind(PETSC_STATIC_ARRAY_LENGTH(PetscCIOptions), PetscCIOptions, name, &idx, &found));
922660278c0SBarry Smith   return found;
923660278c0SBarry Smith }
924660278c0SBarry Smith 
925e5c89e4eSSatish Balay /*@C
92688c29154SBarry Smith   PetscOptionsView - Prints the options that have been loaded. This is
927e5c89e4eSSatish Balay   useful for debugging purposes.
928e5c89e4eSSatish Balay 
929c3339decSBarry Smith   Logically Collective
930e5c89e4eSSatish Balay 
931d8d19677SJose E. Roman   Input Parameters:
93220f4b53cSBarry Smith + options - options database, use `NULL` for default global database
933811af0c4SBarry Smith - viewer  - must be an `PETSCVIEWERASCII` viewer
934e5c89e4eSSatish Balay 
935e5c89e4eSSatish Balay   Options Database Key:
936811af0c4SBarry Smith . -options_view - Activates `PetscOptionsView()` within `PetscFinalize()`
937e5c89e4eSSatish Balay 
93820f4b53cSBarry Smith   Level: advanced
93920f4b53cSBarry Smith 
940811af0c4SBarry Smith   Note:
94121532e8aSBarry Smith   Only the MPI rank 0 of the `MPI_Comm` used to create view prints the option values. Other processes
9421c9f3c13SBarry Smith   may have different values but they are not printed.
9431c9f3c13SBarry Smith 
944db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`
945e5c89e4eSSatish Balay @*/
946d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsView(PetscOptions options, PetscViewer viewer)
947d71ae5a4SJacob Faibussowitsch {
948660278c0SBarry Smith   PetscInt  i, N = 0;
94988c29154SBarry Smith   PetscBool isascii;
950e5c89e4eSSatish Balay 
951e5c89e4eSSatish Balay   PetscFunctionBegin;
9522d747510SLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
953c5929fdfSBarry Smith   options = options ? options : defaultoptions;
95488c29154SBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
9559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
95628b400f6SJacob Faibussowitsch   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Only supports ASCII viewer");
95788c29154SBarry Smith 
958660278c0SBarry Smith   for (i = 0; i < options->N; i++) {
959660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
960660278c0SBarry Smith     N++;
961660278c0SBarry Smith   }
962660278c0SBarry Smith 
963660278c0SBarry Smith   if (!N) {
9649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "#No PETSc Option Table entries\n"));
9653ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
96630694fe9SBarry Smith   }
9672d747510SLisandro Dalcin 
9689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "#PETSc Option Table entries:\n"));
969e5c89e4eSSatish Balay   for (i = 0; i < options->N; i++) {
970660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
971e5c89e4eSSatish Balay     if (options->values[i]) {
9729355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s %s", options->names[i], options->values[i]));
973e5c89e4eSSatish Balay     } else {
9749355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s", options->names[i]));
975e5c89e4eSSatish Balay     }
9769355ec05SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, " # (source: %s)\n", PetscOptionSources[options->source[i]]));
977e5c89e4eSSatish Balay   }
9789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "#End of PETSc Option Table entries\n"));
9793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
980e5c89e4eSSatish Balay }
981e5c89e4eSSatish Balay 
982e11779c2SBarry Smith /*
983e11779c2SBarry Smith    Called by error handlers to print options used in run
984e11779c2SBarry Smith */
985d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftError(void)
986d71ae5a4SJacob Faibussowitsch {
987f4bc716fSBarry Smith   PetscInt i, nopt = 0;
988f4bc716fSBarry Smith 
989f4bc716fSBarry Smith   for (i = 0; i < defaultoptions->N; i++) {
990f4bc716fSBarry Smith     if (!defaultoptions->used[i]) {
991f4bc716fSBarry Smith       if (PetscCIOption(defaultoptions->names[i])) continue;
992f4bc716fSBarry Smith       nopt++;
993f4bc716fSBarry Smith     }
994f4bc716fSBarry Smith   }
995f4bc716fSBarry Smith   if (nopt) {
9967d44c3adSBarry Smith     PetscCall((*PetscErrorPrintf)("WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!\n"));
997f4bc716fSBarry Smith     for (i = 0; i < defaultoptions->N; i++) {
998f4bc716fSBarry Smith       if (!defaultoptions->used[i]) {
999f4bc716fSBarry Smith         if (PetscCIOption(defaultoptions->names[i])) continue;
10003ba16761SJacob 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]]));
10013ba16761SJacob Faibussowitsch         else PetscCall((*PetscErrorPrintf)("  Option left: name:-%s (no value) source: %s\n", defaultoptions->names[i], PetscOptionSources[defaultoptions->source[i]]));
1002f4bc716fSBarry Smith       }
1003f4bc716fSBarry Smith     }
1004f4bc716fSBarry Smith   }
10053ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1006f4bc716fSBarry Smith }
1007f4bc716fSBarry Smith 
1008d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
1009d71ae5a4SJacob Faibussowitsch {
1010660278c0SBarry Smith   PetscInt     i, N = 0;
10114416b707SBarry Smith   PetscOptions options = defaultoptions;
1012e11779c2SBarry Smith 
1013660278c0SBarry Smith   for (i = 0; i < options->N; i++) {
1014660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
1015660278c0SBarry Smith     N++;
1016660278c0SBarry Smith   }
1017660278c0SBarry Smith 
1018660278c0SBarry Smith   if (N) {
10193ba16761SJacob Faibussowitsch     PetscCall((*PetscErrorPrintf)("PETSc Option Table entries:\n"));
1020e11779c2SBarry Smith   } else {
10213ba16761SJacob Faibussowitsch     PetscCall((*PetscErrorPrintf)("No PETSc Option Table entries\n"));
1022e11779c2SBarry Smith   }
1023e11779c2SBarry Smith   for (i = 0; i < options->N; i++) {
1024660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
1025e11779c2SBarry Smith     if (options->values[i]) {
10263ba16761SJacob Faibussowitsch       PetscCall((*PetscErrorPrintf)("-%s %s (source: %s)\n", options->names[i], options->values[i], PetscOptionSources[options->source[i]]));
1027e11779c2SBarry Smith     } else {
10283ba16761SJacob Faibussowitsch       PetscCall((*PetscErrorPrintf)("-%s (source: %s)\n", options->names[i], PetscOptionSources[options->source[i]]));
1029e11779c2SBarry Smith     }
1030e11779c2SBarry Smith   }
10313ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
1032e11779c2SBarry Smith }
1033e11779c2SBarry Smith 
1034e5c89e4eSSatish Balay /*@C
103574e0666dSJed Brown   PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.
103674e0666dSJed Brown 
10371c9f3c13SBarry Smith   Logically Collective
103874e0666dSJed Brown 
1039d8d19677SJose E. Roman   Input Parameters:
104020f4b53cSBarry Smith + options - options database, or `NULL` for the default global database
1041c5929fdfSBarry Smith - prefix  - The string to append to the existing prefix
10429db968c8SJed Brown 
10439db968c8SJed Brown   Options Database Keys:
10449db968c8SJed Brown + -prefix_push <some_prefix_> - push the given prefix
10459db968c8SJed Brown - -prefix_pop                 - pop the last prefix
10469db968c8SJed Brown 
104720f4b53cSBarry Smith   Level: advanced
104820f4b53cSBarry Smith 
10499db968c8SJed Brown   Notes:
105021532e8aSBarry Smith   It is common to use this in conjunction with `-options_file` as in
10519314d9b7SBarry Smith .vb
10529314d9b7SBarry Smith  -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
10539314d9b7SBarry Smith .ve
105421532e8aSBarry Smith   where the files no longer require all options to be prefixed with `-system2_`.
105574e0666dSJed Brown 
105621532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
10571c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
10581c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
10591c9f3c13SBarry Smith   on different ranks.
10601c9f3c13SBarry Smith 
1061db781477SPatrick Sanan .seealso: `PetscOptionsPrefixPop()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
106274e0666dSJed Brown @*/
1063d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPrefixPush(PetscOptions options, const char prefix[])
1064d71ae5a4SJacob Faibussowitsch {
106574e0666dSJed Brown   size_t    n;
106674e0666dSJed Brown   PetscInt  start;
10679355ec05SMatthew G. Knepley   char      key[PETSC_MAX_OPTION_NAME + 1];
10682d747510SLisandro Dalcin   PetscBool valid;
106974e0666dSJed Brown 
107074e0666dSJed Brown   PetscFunctionBegin;
10714f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
1072c5929fdfSBarry Smith   options = options ? options : defaultoptions;
107300045ab3SPierre Jolivet   PetscCheck(options->prefixind < MAXPREFIXES, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum depth of prefix stack %d exceeded, recompile src/sys/objects/options.c with larger value for MAXPREFIXES", MAXPREFIXES);
10742d747510SLisandro Dalcin   key[0] = '-'; /* keys must start with '-' */
10759566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(key + 1, prefix, sizeof(key) - 1));
10769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(key, &valid));
10778bf569ecSLisandro 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 */
107828b400f6SJacob 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" : "");
107974e0666dSJed Brown   start = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
10809566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(prefix, &n));
108108401ef6SPierre Jolivet   PetscCheck(n + 1 <= sizeof(options->prefix) - start, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum prefix length %zu exceeded", sizeof(options->prefix));
10829566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(options->prefix + start, prefix, n + 1));
108374e0666dSJed Brown   options->prefixstack[options->prefixind++] = start + n;
10843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108574e0666dSJed Brown }
108674e0666dSJed Brown 
1087c5929fdfSBarry Smith /*@C
1088811af0c4SBarry Smith   PetscOptionsPrefixPop - Remove the latest options prefix, see `PetscOptionsPrefixPush()` for details
108974e0666dSJed Brown 
1090811af0c4SBarry Smith   Logically Collective on the `MPI_Comm` used when called `PetscOptionsPrefixPush()`
109174e0666dSJed Brown 
1092811af0c4SBarry Smith   Input Parameter:
109320f4b53cSBarry Smith . options - options database, or `NULL` for the default global database
1094c5929fdfSBarry Smith 
109574e0666dSJed Brown   Level: advanced
109674e0666dSJed Brown 
1097db781477SPatrick Sanan .seealso: `PetscOptionsPrefixPush()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
109874e0666dSJed Brown @*/
1099d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1100d71ae5a4SJacob Faibussowitsch {
110174e0666dSJed Brown   PetscInt offset;
110274e0666dSJed Brown 
110374e0666dSJed Brown   PetscFunctionBegin;
1104c5929fdfSBarry Smith   options = options ? options : defaultoptions;
110508401ef6SPierre Jolivet   PetscCheck(options->prefixind >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "More prefixes popped than pushed");
110674e0666dSJed Brown   options->prefixind--;
110774e0666dSJed Brown   offset                  = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
110874e0666dSJed Brown   options->prefix[offset] = 0;
11093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111074e0666dSJed Brown }
111174e0666dSJed Brown 
1112a542b6e8SBarry Smith /*@C
1113a542b6e8SBarry Smith   PetscOptionsClear - Removes all options form the database leaving it empty.
1114a542b6e8SBarry Smith 
11151c9f3c13SBarry Smith   Logically Collective
11161c9f3c13SBarry Smith 
1117811af0c4SBarry Smith   Input Parameter:
111820f4b53cSBarry Smith . options - options database, use `NULL` for the default global database
1119c5929fdfSBarry Smith 
112020f4b53cSBarry Smith   Level: developer
112120f4b53cSBarry Smith 
112220f4b53cSBarry Smith   Note:
112321532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
11241c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
11251c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
11261c9f3c13SBarry Smith   on different ranks.
11271c9f3c13SBarry Smith 
1128db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`
1129a542b6e8SBarry Smith @*/
1130d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsClear(PetscOptions options)
1131d71ae5a4SJacob Faibussowitsch {
1132a542b6e8SBarry Smith   PetscInt i;
1133a542b6e8SBarry Smith 
113439a651e2SJacob Faibussowitsch   PetscFunctionBegin;
1135c5929fdfSBarry Smith   options = options ? options : defaultoptions;
11363ba16761SJacob Faibussowitsch   if (!options) PetscFunctionReturn(PETSC_SUCCESS);
11372d747510SLisandro Dalcin 
1138a542b6e8SBarry Smith   for (i = 0; i < options->N; i++) {
1139a542b6e8SBarry Smith     if (options->names[i]) free(options->names[i]);
1140a542b6e8SBarry Smith     if (options->values[i]) free(options->values[i]);
1141a542b6e8SBarry Smith   }
11422d747510SLisandro Dalcin   options->N = 0;
11439355ec05SMatthew G. Knepley   free(options->names);
11449355ec05SMatthew G. Knepley   free(options->values);
11459355ec05SMatthew G. Knepley   free(options->used);
11469355ec05SMatthew G. Knepley   free(options->source);
11479355ec05SMatthew G. Knepley   options->names  = NULL;
11489355ec05SMatthew G. Knepley   options->values = NULL;
11499355ec05SMatthew G. Knepley   options->used   = NULL;
11509355ec05SMatthew G. Knepley   options->source = NULL;
11519355ec05SMatthew G. Knepley   options->Nalloc = 0;
11522d747510SLisandro Dalcin 
11539355ec05SMatthew G. Knepley   for (i = 0; i < options->Na; i++) {
1154a542b6e8SBarry Smith     free(options->aliases1[i]);
1155a542b6e8SBarry Smith     free(options->aliases2[i]);
1156a542b6e8SBarry Smith   }
11579355ec05SMatthew G. Knepley   options->Na = 0;
11589355ec05SMatthew G. Knepley   free(options->aliases1);
11599355ec05SMatthew G. Knepley   free(options->aliases2);
11609355ec05SMatthew G. Knepley   options->aliases1 = options->aliases2 = NULL;
11619355ec05SMatthew G. Knepley   options->Naalloc                      = 0;
1162a542b6e8SBarry Smith 
11632d747510SLisandro Dalcin   /* destroy hash table */
11642d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
11652d747510SLisandro Dalcin   options->ht = NULL;
11660eb63584SBarry Smith 
11672d747510SLisandro Dalcin   options->prefixind  = 0;
11682d747510SLisandro Dalcin   options->prefix[0]  = 0;
11692d747510SLisandro Dalcin   options->help       = PETSC_FALSE;
11709355ec05SMatthew G. Knepley   options->help_intro = PETSC_FALSE;
11713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11724416b707SBarry Smith }
11734416b707SBarry Smith 
11742d747510SLisandro Dalcin /*@C
11752d747510SLisandro Dalcin   PetscOptionsSetAlias - Makes a key and alias for another key
11762d747510SLisandro Dalcin 
11771c9f3c13SBarry Smith   Logically Collective
11782d747510SLisandro Dalcin 
11792d747510SLisandro Dalcin   Input Parameters:
118020f4b53cSBarry Smith + options - options database, or `NULL` for default global database
11812d747510SLisandro Dalcin . newname - the alias
11822d747510SLisandro Dalcin - oldname - the name that alias will refer to
11832d747510SLisandro Dalcin 
11842d747510SLisandro Dalcin   Level: advanced
11852d747510SLisandro Dalcin 
118620f4b53cSBarry Smith   Note:
118721532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
11881c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
11891c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
11901c9f3c13SBarry Smith   on different ranks.
11911c9f3c13SBarry Smith 
1192c2e3fba1SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1193c2e3fba1SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1194db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1195c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1196db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1197db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
11982d747510SLisandro Dalcin @*/
1199d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSetAlias(PetscOptions options, const char newname[], const char oldname[])
1200d71ae5a4SJacob Faibussowitsch {
12012d747510SLisandro Dalcin   size_t    len;
12029210b8eaSVaclav Hapla   PetscBool valid;
12032d747510SLisandro Dalcin 
12042d747510SLisandro Dalcin   PetscFunctionBegin;
12054f572ea9SToby Isaac   PetscAssertPointer(newname, 2);
12064f572ea9SToby Isaac   PetscAssertPointer(oldname, 3);
12072d747510SLisandro Dalcin   options = options ? options : defaultoptions;
12089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(newname, &valid));
120928b400f6SJacob Faibussowitsch   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliased option %s", newname);
12109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsValidKey(oldname, &valid));
121128b400f6SJacob Faibussowitsch   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliasee option %s", oldname);
12122d747510SLisandro Dalcin 
12139355ec05SMatthew G. Knepley   if (options->Na == options->Naalloc) {
12149355ec05SMatthew G. Knepley     char **tmpA1, **tmpA2;
12152d747510SLisandro Dalcin 
12169355ec05SMatthew G. Knepley     options->Naalloc = PetscMax(4, options->Naalloc * 2);
12179355ec05SMatthew G. Knepley     tmpA1            = (char **)malloc(options->Naalloc * sizeof(char *));
12189355ec05SMatthew G. Knepley     tmpA2            = (char **)malloc(options->Naalloc * sizeof(char *));
12199355ec05SMatthew G. Knepley     for (int i = 0; i < options->Na; ++i) {
12209355ec05SMatthew G. Knepley       tmpA1[i] = options->aliases1[i];
12219355ec05SMatthew G. Knepley       tmpA2[i] = options->aliases2[i];
12229355ec05SMatthew G. Knepley     }
12239355ec05SMatthew G. Knepley     free(options->aliases1);
12249355ec05SMatthew G. Knepley     free(options->aliases2);
12259355ec05SMatthew G. Knepley     options->aliases1 = tmpA1;
12269355ec05SMatthew G. Knepley     options->aliases2 = tmpA2;
12279355ec05SMatthew G. Knepley   }
12289371c9d4SSatish Balay   newname++;
12299371c9d4SSatish Balay   oldname++;
12309566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(newname, &len));
12319355ec05SMatthew G. Knepley   options->aliases1[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1232c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(options->aliases1[options->Na], newname, len + 1));
12339566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(oldname, &len));
12349355ec05SMatthew G. Knepley   options->aliases2[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1235c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(options->aliases2[options->Na], oldname, len + 1));
12369355ec05SMatthew G. Knepley   ++options->Na;
12373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12382d747510SLisandro Dalcin }
12394416b707SBarry Smith 
1240e5c89e4eSSatish Balay /*@C
1241e5c89e4eSSatish Balay   PetscOptionsSetValue - Sets an option name-value pair in the options
1242e5c89e4eSSatish Balay   database, overriding whatever is already present.
1243e5c89e4eSSatish Balay 
12441c9f3c13SBarry Smith   Logically Collective
1245e5c89e4eSSatish Balay 
1246e5c89e4eSSatish Balay   Input Parameters:
124720f4b53cSBarry Smith + options - options database, use `NULL` for the default global database
1248c5929fdfSBarry Smith . name    - name of option, this SHOULD have the - prepended
124920f4b53cSBarry Smith - value   - the option value (not used for all options, so can be `NULL`)
1250e5c89e4eSSatish Balay 
1251e5c89e4eSSatish Balay   Level: intermediate
1252e5c89e4eSSatish Balay 
1253e5c89e4eSSatish Balay   Note:
1254811af0c4SBarry Smith   This function can be called BEFORE `PetscInitialize()`
1255d49172ceSBarry Smith 
125621532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
12571c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
12581c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
12591c9f3c13SBarry Smith   on different ranks.
12601c9f3c13SBarry Smith 
1261aec76313SJacob Faibussowitsch   Developer Notes:
126220f4b53cSBarry Smith   Uses malloc() directly because PETSc may not be initialized yet.
1263b0250c70SBarry Smith 
1264db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`, `PetscOptionsClearValue()`
1265e5c89e4eSSatish Balay @*/
1266d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSetValue(PetscOptions options, const char name[], const char value[])
1267d71ae5a4SJacob Faibussowitsch {
126839a651e2SJacob Faibussowitsch   PetscFunctionBegin;
12699355ec05SMatthew G. Knepley   PetscCall(PetscOptionsSetValue_Private(options, name, value, NULL, PETSC_OPT_CODE));
12703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1271c5b5d8d5SVaclav Hapla }
1272c5b5d8d5SVaclav Hapla 
12739355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options, const char name[], const char value[], int *pos, PetscOptionSource source)
1274d71ae5a4SJacob Faibussowitsch {
1275e5c89e4eSSatish Balay   size_t    len;
12769355ec05SMatthew G. Knepley   int       n, i;
1277e5c89e4eSSatish Balay   char    **names;
12789355ec05SMatthew G. Knepley   char      fullname[PETSC_MAX_OPTION_NAME] = "";
1279c5b5d8d5SVaclav Hapla   PetscBool flg;
1280e5c89e4eSSatish Balay 
128139a651e2SJacob Faibussowitsch   PetscFunctionBegin;
12827272c0d2SVaclav Hapla   if (!options) {
12839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsCreateDefault());
12847272c0d2SVaclav Hapla     options = defaultoptions;
1285c5929fdfSBarry Smith   }
128639a651e2SJacob Faibussowitsch   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "name %s must start with '-'", name);
12872d747510SLisandro Dalcin 
12889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsSkipPrecedent(options, name, &flg));
12893ba16761SJacob Faibussowitsch   if (flg) PetscFunctionReturn(PETSC_SUCCESS);
1290e5c89e4eSSatish Balay 
12912d747510SLisandro Dalcin   name++; /* skip starting dash */
12922d747510SLisandro Dalcin 
129374e0666dSJed Brown   if (options->prefixind > 0) {
1294d49172ceSBarry Smith     strncpy(fullname, options->prefix, sizeof(fullname));
12952d747510SLisandro Dalcin     fullname[sizeof(fullname) - 1] = 0;
129689ae1891SBarry Smith     strncat(fullname, name, sizeof(fullname) - strlen(fullname) - 1);
12972d747510SLisandro Dalcin     fullname[sizeof(fullname) - 1] = 0;
129874e0666dSJed Brown     name                           = fullname;
129974e0666dSJed Brown   }
130074e0666dSJed Brown 
130174e0666dSJed Brown   /* check against aliases */
13029355ec05SMatthew G. Knepley   for (i = 0; i < options->Na; i++) {
13032d747510SLisandro Dalcin     int result = PetscOptNameCmp(options->aliases1[i], name);
13049371c9d4SSatish Balay     if (!result) {
13059371c9d4SSatish Balay       name = options->aliases2[i];
13069371c9d4SSatish Balay       break;
13079371c9d4SSatish Balay     }
1308e5c89e4eSSatish Balay   }
1309e5c89e4eSSatish Balay 
13102d747510SLisandro Dalcin   /* slow search */
13119355ec05SMatthew G. Knepley   n     = options->N;
1312e5c89e4eSSatish Balay   names = options->names;
13139355ec05SMatthew G. Knepley   for (i = 0; i < options->N; i++) {
13142d747510SLisandro Dalcin     int result = PetscOptNameCmp(names[i], name);
13152d747510SLisandro Dalcin     if (!result) {
13169371c9d4SSatish Balay       n = i;
13179371c9d4SSatish Balay       goto setvalue;
13182d747510SLisandro Dalcin     } else if (result > 0) {
13199371c9d4SSatish Balay       n = i;
13209371c9d4SSatish Balay       break;
1321e5c89e4eSSatish Balay     }
1322e5c89e4eSSatish Balay   }
13239355ec05SMatthew G. Knepley   if (options->N == options->Nalloc) {
13249355ec05SMatthew G. Knepley     char             **names, **values;
13259355ec05SMatthew G. Knepley     PetscBool         *used;
13269355ec05SMatthew G. Knepley     PetscOptionSource *source;
13279355ec05SMatthew G. Knepley 
13289355ec05SMatthew G. Knepley     options->Nalloc = PetscMax(10, options->Nalloc * 2);
13299355ec05SMatthew G. Knepley     names           = (char **)malloc(options->Nalloc * sizeof(char *));
13309355ec05SMatthew G. Knepley     values          = (char **)malloc(options->Nalloc * sizeof(char *));
13319355ec05SMatthew G. Knepley     used            = (PetscBool *)malloc(options->Nalloc * sizeof(PetscBool));
13329355ec05SMatthew G. Knepley     source          = (PetscOptionSource *)malloc(options->Nalloc * sizeof(PetscOptionSource));
13339355ec05SMatthew G. Knepley     for (int i = 0; i < options->N; ++i) {
13349355ec05SMatthew G. Knepley       names[i]  = options->names[i];
13359355ec05SMatthew G. Knepley       values[i] = options->values[i];
13369355ec05SMatthew G. Knepley       used[i]   = options->used[i];
13379355ec05SMatthew G. Knepley       source[i] = options->source[i];
13389355ec05SMatthew G. Knepley     }
13399355ec05SMatthew G. Knepley     free(options->names);
13409355ec05SMatthew G. Knepley     free(options->values);
13419355ec05SMatthew G. Knepley     free(options->used);
13429355ec05SMatthew G. Knepley     free(options->source);
13439355ec05SMatthew G. Knepley     options->names  = names;
13449355ec05SMatthew G. Knepley     options->values = values;
13459355ec05SMatthew G. Knepley     options->used   = used;
13469355ec05SMatthew G. Knepley     options->source = source;
13479355ec05SMatthew G. Knepley   }
134839a651e2SJacob Faibussowitsch 
13492d747510SLisandro Dalcin   /* shift remaining values up 1 */
13509355ec05SMatthew G. Knepley   for (i = options->N; i > n; i--) {
13515e8c5e88SLisandro Dalcin     options->names[i]  = options->names[i - 1];
1352e5c89e4eSSatish Balay     options->values[i] = options->values[i - 1];
1353e5c89e4eSSatish Balay     options->used[i]   = options->used[i - 1];
13549355ec05SMatthew G. Knepley     options->source[i] = options->source[i - 1];
1355e5c89e4eSSatish Balay   }
13562d747510SLisandro Dalcin   options->names[n]  = NULL;
13572d747510SLisandro Dalcin   options->values[n] = NULL;
13582d747510SLisandro Dalcin   options->used[n]   = PETSC_FALSE;
13599355ec05SMatthew G. Knepley   options->source[n] = PETSC_OPT_CODE;
13602d747510SLisandro Dalcin   options->N++;
13612d747510SLisandro Dalcin 
13622d747510SLisandro Dalcin   /* destroy hash table */
13632d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
13642d747510SLisandro Dalcin   options->ht = NULL;
13652d747510SLisandro Dalcin 
13662d747510SLisandro Dalcin   /* set new name */
136770d8d27cSBarry Smith   len               = strlen(name);
13685e8c5e88SLisandro Dalcin   options->names[n] = (char *)malloc((len + 1) * sizeof(char));
136939a651e2SJacob Faibussowitsch   PetscCheck(options->names[n], PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate option name");
1370d49172ceSBarry Smith   strcpy(options->names[n], name);
13712d747510SLisandro Dalcin 
13722d747510SLisandro Dalcin setvalue:
13732d747510SLisandro Dalcin   /* set new value */
13742d747510SLisandro Dalcin   if (options->values[n]) free(options->values[n]);
1375d49172ceSBarry Smith   len = value ? strlen(value) : 0;
13765e8c5e88SLisandro Dalcin   if (len) {
1377e5c89e4eSSatish Balay     options->values[n] = (char *)malloc((len + 1) * sizeof(char));
1378d49172ceSBarry Smith     if (!options->values[n]) return PETSC_ERR_MEM;
1379d49172ceSBarry Smith     strcpy(options->values[n], value);
1380*432b765aSRené Chenard     options->values[n][len] = '\0';
13812d747510SLisandro Dalcin   } else {
13822d747510SLisandro Dalcin     options->values[n] = NULL;
13832d747510SLisandro Dalcin   }
13849355ec05SMatthew G. Knepley   options->source[n] = source;
13852d747510SLisandro Dalcin 
138691ad3481SVaclav Hapla   /* handle -help so that it can be set from anywhere */
138791ad3481SVaclav Hapla   if (!PetscOptNameCmp(name, "help")) {
138891ad3481SVaclav Hapla     options->help       = PETSC_TRUE;
1389d06005cbSLisandro Dalcin     options->help_intro = (value && !PetscOptNameCmp(value, "intro")) ? PETSC_TRUE : PETSC_FALSE;
139091ad3481SVaclav Hapla     options->used[n]    = PETSC_TRUE;
139191ad3481SVaclav Hapla   }
139291ad3481SVaclav Hapla 
1393*432b765aSRené Chenard   PetscCall(PetscOptionsMonitor(options, name, value ? value : "", source));
1394c5b5d8d5SVaclav Hapla   if (pos) *pos = n;
13953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1396e5c89e4eSSatish Balay }
1397e5c89e4eSSatish Balay 
1398e5c89e4eSSatish Balay /*@C
1399e5c89e4eSSatish Balay   PetscOptionsClearValue - Clears an option name-value pair in the options
1400e5c89e4eSSatish Balay   database, overriding whatever is already present.
1401e5c89e4eSSatish Balay 
14021c9f3c13SBarry Smith   Logically Collective
1403e5c89e4eSSatish Balay 
1404d8d19677SJose E. Roman   Input Parameters:
140520f4b53cSBarry Smith + options - options database, use `NULL` for the default global database
1406a2b725a8SWilliam Gropp - name    - name of option, this SHOULD have the - prepended
1407e5c89e4eSSatish Balay 
1408e5c89e4eSSatish Balay   Level: intermediate
1409e5c89e4eSSatish Balay 
1410811af0c4SBarry Smith   Note:
141121532e8aSBarry Smith   The collectivity of this routine is complex; only the MPI processes that call this routine will
14121c9f3c13SBarry Smith   have the affect of these options. If some processes that create objects call this routine and others do
14131c9f3c13SBarry Smith   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
14141c9f3c13SBarry Smith   on different ranks.
14151c9f3c13SBarry Smith 
1416db781477SPatrick Sanan .seealso: `PetscOptionsInsert()`
1417e5c89e4eSSatish Balay @*/
1418d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsClearValue(PetscOptions options, const char name[])
1419d71ae5a4SJacob Faibussowitsch {
14202d747510SLisandro Dalcin   int    N, n, i;
14212d747510SLisandro Dalcin   char **names;
1422e5c89e4eSSatish Balay 
1423e5c89e4eSSatish Balay   PetscFunctionBegin;
1424c5929fdfSBarry Smith   options = options ? options : defaultoptions;
1425cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1426c9dcd962SLisandro Dalcin   if (!PetscOptNameCmp(name, "-help")) options->help = options->help_intro = PETSC_FALSE;
14272d747510SLisandro Dalcin 
14282d747510SLisandro Dalcin   name++; /* skip starting dash */
14292d747510SLisandro Dalcin 
14302d747510SLisandro Dalcin   /* slow search */
14312d747510SLisandro Dalcin   N = n = options->N;
1432e5c89e4eSSatish Balay   names = options->names;
1433e5c89e4eSSatish Balay   for (i = 0; i < N; i++) {
14342d747510SLisandro Dalcin     int result = PetscOptNameCmp(names[i], name);
14352d747510SLisandro Dalcin     if (!result) {
14369371c9d4SSatish Balay       n = i;
14379371c9d4SSatish Balay       break;
14382d747510SLisandro Dalcin     } else if (result > 0) {
14399371c9d4SSatish Balay       n = N;
14409371c9d4SSatish Balay       break;
1441e5c89e4eSSatish Balay     }
14422d747510SLisandro Dalcin   }
14433ba16761SJacob Faibussowitsch   if (n == N) PetscFunctionReturn(PETSC_SUCCESS); /* it was not present */
1444e5c89e4eSSatish Balay 
14452d747510SLisandro Dalcin   /* remove name and value */
14462d747510SLisandro Dalcin   if (options->names[n]) free(options->names[n]);
14472d747510SLisandro Dalcin   if (options->values[n]) free(options->values[n]);
1448e5c89e4eSSatish Balay   /* shift remaining values down 1 */
1449e5c89e4eSSatish Balay   for (i = n; i < N - 1; i++) {
14505e8c5e88SLisandro Dalcin     options->names[i]  = options->names[i + 1];
1451e5c89e4eSSatish Balay     options->values[i] = options->values[i + 1];
1452e5c89e4eSSatish Balay     options->used[i]   = options->used[i + 1];
14539355ec05SMatthew G. Knepley     options->source[i] = options->source[i + 1];
1454e5c89e4eSSatish Balay   }
1455e5c89e4eSSatish Balay   options->N--;
14562d747510SLisandro Dalcin 
14572d747510SLisandro Dalcin   /* destroy hash table */
14582d747510SLisandro Dalcin   kh_destroy(HO, options->ht);
14592d747510SLisandro Dalcin   options->ht = NULL;
14602d747510SLisandro Dalcin 
14619355ec05SMatthew G. Knepley   PetscCall(PetscOptionsMonitor(options, name, NULL, PETSC_OPT_CODE));
14623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1463e5c89e4eSSatish Balay }
1464e5c89e4eSSatish Balay 
1465e5c89e4eSSatish Balay /*@C
14662d747510SLisandro Dalcin   PetscOptionsFindPair - Gets an option name-value pair from the options database.
1467e5c89e4eSSatish Balay 
14682d747510SLisandro Dalcin   Not Collective
1469e5c89e4eSSatish Balay 
1470e5c89e4eSSatish Balay   Input Parameters:
147120f4b53cSBarry Smith + options - options database, use `NULL` for the default global database
147220f4b53cSBarry Smith . pre     - the string to prepend to the name or `NULL`, this SHOULD NOT have the "-" prepended
14732d747510SLisandro Dalcin - name    - name of option, this SHOULD have the "-" prepended
1474e5c89e4eSSatish Balay 
14752d747510SLisandro Dalcin   Output Parameters:
14762d747510SLisandro Dalcin + value - the option value (optional, not used for all options)
14772d747510SLisandro Dalcin - set   - whether the option is set (optional)
1478e5c89e4eSSatish Balay 
147920f4b53cSBarry Smith   Level: developer
148020f4b53cSBarry Smith 
1481811af0c4SBarry Smith   Note:
14829666a313SBarry Smith   Each process may find different values or no value depending on how options were inserted into the database
14831c9f3c13SBarry Smith 
1484db781477SPatrick Sanan .seealso: `PetscOptionsSetValue()`, `PetscOptionsClearValue()`
1485e5c89e4eSSatish Balay @*/
1486d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsFindPair(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1487d71ae5a4SJacob Faibussowitsch {
14889355ec05SMatthew G. Knepley   char      buf[PETSC_MAX_OPTION_NAME];
1489daabea38SBarry Smith   PetscBool usehashtable = PETSC_TRUE;
14902d747510SLisandro Dalcin   PetscBool matchnumbers = PETSC_TRUE;
1491e5c89e4eSSatish Balay 
1492e5c89e4eSSatish Balay   PetscFunctionBegin;
1493c5929fdfSBarry Smith   options = options ? options : defaultoptions;
149408401ef6SPierre Jolivet   PetscCheck(!pre || !PetscUnlikely(pre[0] == '-'), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1495cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1496e5c89e4eSSatish Balay 
14972d747510SLisandro Dalcin   name++; /* skip starting dash */
1498e5c89e4eSSatish Balay 
14997cd08cecSJed Brown   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
15002d747510SLisandro Dalcin   if (pre && pre[0]) {
15012d747510SLisandro Dalcin     char *ptr = buf;
15029371c9d4SSatish Balay     if (name[0] == '-') {
15039371c9d4SSatish Balay       *ptr++ = '-';
15049371c9d4SSatish Balay       name++;
15059371c9d4SSatish Balay     }
15069566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(ptr, pre, buf + sizeof(buf) - ptr));
15079566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
15082d747510SLisandro Dalcin     name = buf;
15097cd08cecSJed Brown   }
15102d747510SLisandro Dalcin 
151176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
15122f828895SJed Brown     PetscBool valid;
15139355ec05SMatthew G. Knepley     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
15149566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
15159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(key, &valid));
151628b400f6SJacob Faibussowitsch     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
15172f828895SJed Brown   }
1518e5c89e4eSSatish Balay 
15192d747510SLisandro Dalcin   if (!options->ht && usehashtable) {
15202d747510SLisandro Dalcin     int          i, ret;
15212d747510SLisandro Dalcin     khiter_t     it;
15222d747510SLisandro Dalcin     khash_t(HO) *ht;
15232d747510SLisandro Dalcin     ht = kh_init(HO);
152428b400f6SJacob Faibussowitsch     PetscCheck(ht, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15252d747510SLisandro Dalcin     ret = kh_resize(HO, ht, options->N * 2); /* twice the required size to reduce risk of collisions */
152628b400f6SJacob Faibussowitsch     PetscCheck(!ret, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15272d747510SLisandro Dalcin     for (i = 0; i < options->N; i++) {
15282d747510SLisandro Dalcin       it = kh_put(HO, ht, options->names[i], &ret);
152908401ef6SPierre Jolivet       PetscCheck(ret == 1, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
15302d747510SLisandro Dalcin       kh_val(ht, it) = i;
15312d747510SLisandro Dalcin     }
15322d747510SLisandro Dalcin     options->ht = ht;
15332d747510SLisandro Dalcin   }
15342d747510SLisandro Dalcin 
15359371c9d4SSatish Balay   if (usehashtable) { /* fast search */
15362d747510SLisandro Dalcin     khash_t(HO) *ht = options->ht;
15372d747510SLisandro Dalcin     khiter_t     it = kh_get(HO, ht, name);
15382d747510SLisandro Dalcin     if (it != kh_end(ht)) {
15392d747510SLisandro Dalcin       int i            = kh_val(ht, it);
1540e5c89e4eSSatish Balay       options->used[i] = PETSC_TRUE;
15412d747510SLisandro Dalcin       if (value) *value = options->values[i];
15422d747510SLisandro Dalcin       if (set) *set = PETSC_TRUE;
15433ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
15442d747510SLisandro Dalcin     }
15459371c9d4SSatish Balay   } else { /* slow search */
15462d747510SLisandro Dalcin     int i, N = options->N;
15472d747510SLisandro Dalcin     for (i = 0; i < N; i++) {
1548daabea38SBarry Smith       int result = PetscOptNameCmp(options->names[i], name);
15492d747510SLisandro Dalcin       if (!result) {
15502d747510SLisandro Dalcin         options->used[i] = PETSC_TRUE;
15512d747510SLisandro Dalcin         if (value) *value = options->values[i];
15522d747510SLisandro Dalcin         if (set) *set = PETSC_TRUE;
15533ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
15542d747510SLisandro Dalcin       } else if (result > 0) {
1555e5c89e4eSSatish Balay         break;
1556e5c89e4eSSatish Balay       }
1557e5c89e4eSSatish Balay     }
15582d747510SLisandro Dalcin   }
15592d747510SLisandro Dalcin 
15602d747510SLisandro Dalcin   /*
15612d747510SLisandro Dalcin    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
15622d747510SLisandro Dalcin    Maybe this special lookup mode should be enabled on request with a push/pop API.
15632d747510SLisandro Dalcin    The feature of matching _%d_ used sparingly in the codebase.
15642d747510SLisandro Dalcin    */
15652d747510SLisandro Dalcin   if (matchnumbers) {
15662d747510SLisandro Dalcin     int i, j, cnt = 0, locs[16], loce[16];
1567e5c89e4eSSatish Balay     /* determine the location and number of all _%d_ in the key */
15682d747510SLisandro Dalcin     for (i = 0; name[i]; i++) {
15692d747510SLisandro Dalcin       if (name[i] == '_') {
15702d747510SLisandro Dalcin         for (j = i + 1; name[j]; j++) {
15712d747510SLisandro Dalcin           if (name[j] >= '0' && name[j] <= '9') continue;
15722d747510SLisandro Dalcin           if (name[j] == '_' && j > i + 1) { /* found a number */
1573e5c89e4eSSatish Balay             locs[cnt]   = i + 1;
1574e5c89e4eSSatish Balay             loce[cnt++] = j + 1;
1575e5c89e4eSSatish Balay           }
15762d747510SLisandro Dalcin           i = j - 1;
1577e5c89e4eSSatish Balay           break;
1578e5c89e4eSSatish Balay         }
1579e5c89e4eSSatish Balay       }
1580e5c89e4eSSatish Balay     }
1581e5c89e4eSSatish Balay     for (i = 0; i < cnt; i++) {
15822d747510SLisandro Dalcin       PetscBool found;
15839355ec05SMatthew G. Knepley       char      opt[PETSC_MAX_OPTION_NAME + 1] = "-", tmp[PETSC_MAX_OPTION_NAME];
15849566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(tmp, name, PetscMin((size_t)(locs[i] + 1), sizeof(tmp))));
15859566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, tmp, sizeof(opt)));
15869566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, name + loce[i], sizeof(opt)));
15879566063dSJacob Faibussowitsch       PetscCall(PetscOptionsFindPair(options, NULL, opt, value, &found));
15889371c9d4SSatish Balay       if (found) {
15899371c9d4SSatish Balay         if (set) *set = PETSC_TRUE;
15903ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
15919371c9d4SSatish Balay       }
1592e5c89e4eSSatish Balay     }
1593e5c89e4eSSatish Balay   }
15942d747510SLisandro Dalcin 
15952d747510SLisandro Dalcin   if (set) *set = PETSC_FALSE;
15963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1597e5c89e4eSSatish Balay }
1598e5c89e4eSSatish Balay 
1599d6ced9c0SMatthew G. Knepley /* Check whether any option begins with pre+name */
160054a546c1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options, const char pre[], const char name[], const char *option[], const char *value[], PetscBool *set)
1601d71ae5a4SJacob Faibussowitsch {
16029355ec05SMatthew G. Knepley   char buf[PETSC_MAX_OPTION_NAME];
1603d6ced9c0SMatthew G. Knepley   int  numCnt = 0, locs[16], loce[16];
1604514bf10dSMatthew G Knepley 
1605514bf10dSMatthew G Knepley   PetscFunctionBegin;
1606c5929fdfSBarry Smith   options = options ? options : defaultoptions;
1607cc73adaaSBarry Smith   PetscCheck(!pre || pre[0] != '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1608cc73adaaSBarry Smith   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1609514bf10dSMatthew G Knepley 
16102d747510SLisandro Dalcin   name++; /* skip starting dash */
1611514bf10dSMatthew G Knepley 
1612514bf10dSMatthew G Knepley   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
16132d747510SLisandro Dalcin   if (pre && pre[0]) {
16142d747510SLisandro Dalcin     char *ptr = buf;
16159371c9d4SSatish Balay     if (name[0] == '-') {
16169371c9d4SSatish Balay       *ptr++ = '-';
16179371c9d4SSatish Balay       name++;
16189371c9d4SSatish Balay     }
16199b15cf9aSJacob Faibussowitsch     PetscCall(PetscStrncpy(ptr, pre, sizeof(buf) - ((ptr == buf) ? 0 : 1)));
16209566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
16212d747510SLisandro Dalcin     name = buf;
1622514bf10dSMatthew G Knepley   }
16232d747510SLisandro Dalcin 
162476bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
1625514bf10dSMatthew G Knepley     PetscBool valid;
16269355ec05SMatthew G. Knepley     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
16279566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
16289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsValidKey(key, &valid));
162928b400f6SJacob Faibussowitsch     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1630514bf10dSMatthew G Knepley   }
1631514bf10dSMatthew G Knepley 
1632d6ced9c0SMatthew G. Knepley   /* determine the location and number of all _%d_ in the key */
1633d6ced9c0SMatthew G. Knepley   {
1634d6ced9c0SMatthew G. Knepley     int i, j;
1635d6ced9c0SMatthew G. Knepley     for (i = 0; name[i]; i++) {
1636d6ced9c0SMatthew G. Knepley       if (name[i] == '_') {
1637d6ced9c0SMatthew G. Knepley         for (j = i + 1; name[j]; j++) {
1638d6ced9c0SMatthew G. Knepley           if (name[j] >= '0' && name[j] <= '9') continue;
1639d6ced9c0SMatthew G. Knepley           if (name[j] == '_' && j > i + 1) { /* found a number */
1640d6ced9c0SMatthew G. Knepley             locs[numCnt]   = i + 1;
1641d6ced9c0SMatthew G. Knepley             loce[numCnt++] = j + 1;
1642d6ced9c0SMatthew G. Knepley           }
1643d6ced9c0SMatthew G. Knepley           i = j - 1;
1644d6ced9c0SMatthew G. Knepley           break;
1645d6ced9c0SMatthew G. Knepley         }
1646d6ced9c0SMatthew G. Knepley       }
1647d6ced9c0SMatthew G. Knepley     }
1648d6ced9c0SMatthew G. Knepley   }
1649d6ced9c0SMatthew G. Knepley 
1650363da2dcSJacob Faibussowitsch   /* slow search */
1651363da2dcSJacob Faibussowitsch   for (int c = -1; c < numCnt; ++c) {
1652363da2dcSJacob Faibussowitsch     char   opt[PETSC_MAX_OPTION_NAME + 2] = "";
16532d747510SLisandro Dalcin     size_t len;
1654d6ced9c0SMatthew G. Knepley 
1655d6ced9c0SMatthew G. Knepley     if (c < 0) {
1656c6a7a370SJeremy L Thompson       PetscCall(PetscStrncpy(opt, name, sizeof(opt)));
1657d6ced9c0SMatthew G. Knepley     } else {
1658363da2dcSJacob Faibussowitsch       PetscCall(PetscStrncpy(opt, name, PetscMin((size_t)(locs[c] + 1), sizeof(opt))));
1659363da2dcSJacob Faibussowitsch       PetscCall(PetscStrlcat(opt, name + loce[c], sizeof(opt) - 1));
1660d6ced9c0SMatthew G. Knepley     }
16619566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(opt, &len));
1662363da2dcSJacob Faibussowitsch     for (int i = 0; i < options->N; i++) {
1663363da2dcSJacob Faibussowitsch       PetscBool match;
1664363da2dcSJacob Faibussowitsch 
16659566063dSJacob Faibussowitsch       PetscCall(PetscStrncmp(options->names[i], opt, len, &match));
1666514bf10dSMatthew G Knepley       if (match) {
1667514bf10dSMatthew G Knepley         options->used[i] = PETSC_TRUE;
166854a546c1SMatthew G. Knepley         if (option) *option = options->names[i];
16692d747510SLisandro Dalcin         if (value) *value = options->values[i];
16702d747510SLisandro Dalcin         if (set) *set = PETSC_TRUE;
16713ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
1672514bf10dSMatthew G Knepley       }
1673514bf10dSMatthew G Knepley     }
16742d747510SLisandro Dalcin   }
16752d747510SLisandro Dalcin 
16762d747510SLisandro Dalcin   if (set) *set = PETSC_FALSE;
16773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1678514bf10dSMatthew G Knepley }
1679514bf10dSMatthew G Knepley 
1680e5c89e4eSSatish Balay /*@C
1681e5c89e4eSSatish Balay   PetscOptionsReject - Generates an error if a certain option is given.
1682e5c89e4eSSatish Balay 
16831c9f3c13SBarry Smith   Not Collective
1684e5c89e4eSSatish Balay 
1685e5c89e4eSSatish Balay   Input Parameters:
168620f4b53cSBarry Smith + options - options database, use `NULL` for default global database
168720f4b53cSBarry Smith . pre     - the option prefix (may be `NULL`)
16882d747510SLisandro Dalcin . name    - the option name one is seeking
168920f4b53cSBarry Smith - mess    - error message (may be `NULL`)
1690e5c89e4eSSatish Balay 
1691e5c89e4eSSatish Balay   Level: advanced
1692e5c89e4eSSatish Balay 
1693c2e3fba1SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1694db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1695db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1696c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1697db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1698db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1699e5c89e4eSSatish Balay @*/
1700d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsReject(PetscOptions options, const char pre[], const char name[], const char mess[])
1701d71ae5a4SJacob Faibussowitsch {
1702ace3abfcSBarry Smith   PetscBool flag = PETSC_FALSE;
1703e5c89e4eSSatish Balay 
1704e5c89e4eSSatish Balay   PetscFunctionBegin;
17059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(options, pre, name, &flag));
1706e5c89e4eSSatish Balay   if (flag) {
170708401ef6SPierre 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);
1708f7d195e4SLawrence Mitchell     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s", pre ? pre : "", name + 1);
1709e5c89e4eSSatish Balay   }
17103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1711e5c89e4eSSatish Balay }
1712e5c89e4eSSatish Balay 
1713e5c89e4eSSatish Balay /*@C
17142d747510SLisandro Dalcin   PetscOptionsHasHelp - Determines whether the "-help" option is in the database.
17152d747510SLisandro Dalcin 
17162d747510SLisandro Dalcin   Not Collective
17172d747510SLisandro Dalcin 
1718811af0c4SBarry Smith   Input Parameter:
171920f4b53cSBarry Smith . options - options database, use `NULL` for default global database
17202d747510SLisandro Dalcin 
1721811af0c4SBarry Smith   Output Parameter:
1722811af0c4SBarry Smith . set - `PETSC_TRUE` if found else `PETSC_FALSE`.
17232d747510SLisandro Dalcin 
17242d747510SLisandro Dalcin   Level: advanced
17252d747510SLisandro Dalcin 
1726db781477SPatrick Sanan .seealso: `PetscOptionsHasName()`
17272d747510SLisandro Dalcin @*/
1728d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasHelp(PetscOptions options, PetscBool *set)
1729d71ae5a4SJacob Faibussowitsch {
17302d747510SLisandro Dalcin   PetscFunctionBegin;
17314f572ea9SToby Isaac   PetscAssertPointer(set, 2);
17322d747510SLisandro Dalcin   options = options ? options : defaultoptions;
17332d747510SLisandro Dalcin   *set    = options->help;
17343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17352d747510SLisandro Dalcin }
17362d747510SLisandro Dalcin 
1737d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options, PetscBool *set)
1738d71ae5a4SJacob Faibussowitsch {
1739d314f959SVaclav Hapla   PetscFunctionBegin;
17404f572ea9SToby Isaac   PetscAssertPointer(set, 2);
1741d314f959SVaclav Hapla   options = options ? options : defaultoptions;
1742d314f959SVaclav Hapla   *set    = options->help_intro;
17433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1744d314f959SVaclav Hapla }
1745d314f959SVaclav Hapla 
17462d747510SLisandro Dalcin /*@C
1747e24fcbf7SPierre 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
1748e24fcbf7SPierre Jolivet   if its value is set to false.
1749e5c89e4eSSatish Balay 
1750e5c89e4eSSatish Balay   Not Collective
1751e5c89e4eSSatish Balay 
1752e5c89e4eSSatish Balay   Input Parameters:
175320f4b53cSBarry Smith + options - options database, use `NULL` for default global database
175420f4b53cSBarry Smith . pre     - string to prepend to the name or `NULL`
17553de71b31SHong Zhang - name    - the option one is seeking
1756e5c89e4eSSatish Balay 
1757811af0c4SBarry Smith   Output Parameter:
1758811af0c4SBarry Smith . set - `PETSC_TRUE` if found else `PETSC_FALSE`.
1759e5c89e4eSSatish Balay 
1760e5c89e4eSSatish Balay   Level: beginner
1761e5c89e4eSSatish Balay 
1762811af0c4SBarry Smith   Note:
1763811af0c4SBarry Smith   In many cases you probably want to use `PetscOptionsGetBool()` instead of calling this, to allowing toggling values.
176490d69ab7SBarry Smith 
1765db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
1766db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1767db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1768c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1769db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1770db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1771e5c89e4eSSatish Balay @*/
1772d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsHasName(PetscOptions options, const char pre[], const char name[], PetscBool *set)
1773d71ae5a4SJacob Faibussowitsch {
17742d747510SLisandro Dalcin   const char *value;
1775ace3abfcSBarry Smith   PetscBool   flag;
1776e5c89e4eSSatish Balay 
1777e5c89e4eSSatish Balay   PetscFunctionBegin;
17789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
177996ef3cdfSSatish Balay   if (set) *set = flag;
17803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1781e5c89e4eSSatish Balay }
1782e5c89e4eSSatish Balay 
1783e5c89e4eSSatish Balay /*@C
17842d747510SLisandro Dalcin   PetscOptionsGetAll - Lists all the options the program was run with in a single string.
17852d747510SLisandro Dalcin 
17862d747510SLisandro Dalcin   Not Collective
17872d747510SLisandro Dalcin 
1788fd292e60Sprj-   Input Parameter:
178920f4b53cSBarry Smith . options - the options database, use `NULL` for the default global database
17902d747510SLisandro Dalcin 
17912d747510SLisandro Dalcin   Output Parameter:
17922d747510SLisandro Dalcin . copts - pointer where string pointer is stored
17932d747510SLisandro Dalcin 
179420f4b53cSBarry Smith   Level: advanced
179520f4b53cSBarry Smith 
17962d747510SLisandro Dalcin   Notes:
1797811af0c4SBarry Smith   The array and each entry in the array should be freed with `PetscFree()`
1798811af0c4SBarry Smith 
17991c9f3c13SBarry Smith   Each process may have different values depending on how the options were inserted into the database
18002d747510SLisandro Dalcin 
1801db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsView()`, `PetscOptionsPush()`, `PetscOptionsPop()`
18022d747510SLisandro Dalcin @*/
1803d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetAll(PetscOptions options, char *copts[])
1804d71ae5a4SJacob Faibussowitsch {
18052d747510SLisandro Dalcin   PetscInt i;
18062d747510SLisandro Dalcin   size_t   len = 1, lent = 0;
18072d747510SLisandro Dalcin   char    *coptions = NULL;
18082d747510SLisandro Dalcin 
18092d747510SLisandro Dalcin   PetscFunctionBegin;
18104f572ea9SToby Isaac   PetscAssertPointer(copts, 2);
18112d747510SLisandro Dalcin   options = options ? options : defaultoptions;
18122d747510SLisandro Dalcin   /* count the length of the required string */
18132d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
18149566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(options->names[i], &lent));
18152d747510SLisandro Dalcin     len += 2 + lent;
18162d747510SLisandro Dalcin     if (options->values[i]) {
18179566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(options->values[i], &lent));
18182d747510SLisandro Dalcin       len += 1 + lent;
18192d747510SLisandro Dalcin     }
18202d747510SLisandro Dalcin   }
18219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &coptions));
18222d747510SLisandro Dalcin   coptions[0] = 0;
18232d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
1824c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, "-", len));
1825c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, options->names[i], len));
1826c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(coptions, " ", len));
18272d747510SLisandro Dalcin     if (options->values[i]) {
1828c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(coptions, options->values[i], len));
1829c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(coptions, " ", len));
18302d747510SLisandro Dalcin     }
18312d747510SLisandro Dalcin   }
18322d747510SLisandro Dalcin   *copts = coptions;
18333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18342d747510SLisandro Dalcin }
18352d747510SLisandro Dalcin 
18362d747510SLisandro Dalcin /*@C
18372d747510SLisandro Dalcin   PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database
18382d747510SLisandro Dalcin 
18392d747510SLisandro Dalcin   Not Collective
18402d747510SLisandro Dalcin 
1841d8d19677SJose E. Roman   Input Parameters:
184220f4b53cSBarry Smith + options - options database, use `NULL` for default global database
18432d747510SLisandro Dalcin - name    - string name of option
18442d747510SLisandro Dalcin 
18452d747510SLisandro Dalcin   Output Parameter:
1846811af0c4SBarry Smith . used - `PETSC_TRUE` if the option was used, otherwise false, including if option was not found in options database
18472d747510SLisandro Dalcin 
18482d747510SLisandro Dalcin   Level: advanced
18492d747510SLisandro Dalcin 
1850811af0c4SBarry Smith   Note:
18519666a313SBarry Smith   The value returned may be different on each process and depends on which options have been processed
18521c9f3c13SBarry Smith   on the given process
18531c9f3c13SBarry Smith 
1854db781477SPatrick Sanan .seealso: `PetscOptionsView()`, `PetscOptionsLeft()`, `PetscOptionsAllUsed()`
18552d747510SLisandro Dalcin @*/
1856d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsUsed(PetscOptions options, const char *name, PetscBool *used)
1857d71ae5a4SJacob Faibussowitsch {
18582d747510SLisandro Dalcin   PetscInt i;
18592d747510SLisandro Dalcin 
18602d747510SLisandro Dalcin   PetscFunctionBegin;
18614f572ea9SToby Isaac   PetscAssertPointer(name, 2);
18624f572ea9SToby Isaac   PetscAssertPointer(used, 3);
18632d747510SLisandro Dalcin   options = options ? options : defaultoptions;
18642d747510SLisandro Dalcin   *used   = PETSC_FALSE;
18652d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
18669566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp(options->names[i], name, used));
18672d747510SLisandro Dalcin     if (*used) {
18682d747510SLisandro Dalcin       *used = options->used[i];
18692d747510SLisandro Dalcin       break;
18702d747510SLisandro Dalcin     }
18712d747510SLisandro Dalcin   }
18723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18732d747510SLisandro Dalcin }
18742d747510SLisandro Dalcin 
1875487a658cSBarry Smith /*@
18762d747510SLisandro Dalcin   PetscOptionsAllUsed - Returns a count of the number of options in the
18772d747510SLisandro Dalcin   database that have never been selected.
18782d747510SLisandro Dalcin 
18792d747510SLisandro Dalcin   Not Collective
18802d747510SLisandro Dalcin 
18812d747510SLisandro Dalcin   Input Parameter:
188220f4b53cSBarry Smith . options - options database, use `NULL` for default global database
18832d747510SLisandro Dalcin 
18842d747510SLisandro Dalcin   Output Parameter:
18852d747510SLisandro Dalcin . N - count of options not used
18862d747510SLisandro Dalcin 
18872d747510SLisandro Dalcin   Level: advanced
18882d747510SLisandro Dalcin 
1889811af0c4SBarry Smith   Note:
18909666a313SBarry Smith   The value returned may be different on each process and depends on which options have been processed
18911c9f3c13SBarry Smith   on the given process
18921c9f3c13SBarry Smith 
1893db781477SPatrick Sanan .seealso: `PetscOptionsView()`
18942d747510SLisandro Dalcin @*/
1895d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsAllUsed(PetscOptions options, PetscInt *N)
1896d71ae5a4SJacob Faibussowitsch {
18972d747510SLisandro Dalcin   PetscInt i, n = 0;
18982d747510SLisandro Dalcin 
18992d747510SLisandro Dalcin   PetscFunctionBegin;
19004f572ea9SToby Isaac   PetscAssertPointer(N, 2);
19012d747510SLisandro Dalcin   options = options ? options : defaultoptions;
19022d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
19032d747510SLisandro Dalcin     if (!options->used[i]) n++;
19042d747510SLisandro Dalcin   }
19052d747510SLisandro Dalcin   *N = n;
19063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19072d747510SLisandro Dalcin }
19082d747510SLisandro Dalcin 
1909487a658cSBarry Smith /*@
19102d747510SLisandro Dalcin   PetscOptionsLeft - Prints to screen any options that were set and never used.
19112d747510SLisandro Dalcin 
19122d747510SLisandro Dalcin   Not Collective
19132d747510SLisandro Dalcin 
19142d747510SLisandro Dalcin   Input Parameter:
191520f4b53cSBarry Smith . options - options database; use `NULL` for default global database
19162d747510SLisandro Dalcin 
19172d747510SLisandro Dalcin   Options Database Key:
1918811af0c4SBarry Smith . -options_left - activates `PetscOptionsAllUsed()` within `PetscFinalize()`
19192d747510SLisandro Dalcin 
192020f4b53cSBarry Smith   Level: advanced
192120f4b53cSBarry Smith 
19223de2bfdfSBarry Smith   Notes:
1923811af0c4SBarry Smith   This is rarely used directly, it is called by `PetscFinalize()` in debug more or if -options_left
19241c9f3c13SBarry Smith   is passed otherwise to help users determine possible mistakes in their usage of options. This
1925811af0c4SBarry Smith   only prints values on process zero of `PETSC_COMM_WORLD`.
1926811af0c4SBarry Smith 
1927811af0c4SBarry Smith   Other processes depending the objects
19281c9f3c13SBarry Smith   used may have different options that are left unused.
19293de2bfdfSBarry Smith 
1930db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`
19312d747510SLisandro Dalcin @*/
1932d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeft(PetscOptions options)
1933d71ae5a4SJacob Faibussowitsch {
19342d747510SLisandro Dalcin   PetscInt     i;
19353de2bfdfSBarry Smith   PetscInt     cnt = 0;
19363de2bfdfSBarry Smith   PetscOptions toptions;
19372d747510SLisandro Dalcin 
19382d747510SLisandro Dalcin   PetscFunctionBegin;
19393de2bfdfSBarry Smith   toptions = options ? options : defaultoptions;
19403de2bfdfSBarry Smith   for (i = 0; i < toptions->N; i++) {
19413de2bfdfSBarry Smith     if (!toptions->used[i]) {
1942660278c0SBarry Smith       if (PetscCIOption(toptions->names[i])) continue;
19433de2bfdfSBarry Smith       if (toptions->values[i]) {
19449355ec05SMatthew 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]]));
19452d747510SLisandro Dalcin       } else {
19469355ec05SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s (no value) source: %s\n", toptions->names[i], PetscOptionSources[toptions->source[i]]));
19472d747510SLisandro Dalcin       }
19482d747510SLisandro Dalcin     }
19492d747510SLisandro Dalcin   }
19503de2bfdfSBarry Smith   if (!options) {
19513de2bfdfSBarry Smith     toptions = defaultoptions;
19523de2bfdfSBarry Smith     while (toptions->previous) {
19533de2bfdfSBarry Smith       cnt++;
19543de2bfdfSBarry Smith       toptions = toptions->previous;
19553de2bfdfSBarry Smith     }
195648a46eb9SPierre 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));
19573de2bfdfSBarry Smith   }
19583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19592d747510SLisandro Dalcin }
19602d747510SLisandro Dalcin 
19612d747510SLisandro Dalcin /*@C
19622d747510SLisandro Dalcin   PetscOptionsLeftGet - Returns all options that were set and never used.
19632d747510SLisandro Dalcin 
19642d747510SLisandro Dalcin   Not Collective
19652d747510SLisandro Dalcin 
19662d747510SLisandro Dalcin   Input Parameter:
196720f4b53cSBarry Smith . options - options database, use `NULL` for default global database
19682d747510SLisandro Dalcin 
1969d8d19677SJose E. Roman   Output Parameters:
1970a2b725a8SWilliam Gropp + N      - count of options not used
19712d747510SLisandro Dalcin . names  - names of options not used
1972a2b725a8SWilliam Gropp - values - values of options not used
19732d747510SLisandro Dalcin 
19742d747510SLisandro Dalcin   Level: advanced
19752d747510SLisandro Dalcin 
19762d747510SLisandro Dalcin   Notes:
1977811af0c4SBarry Smith   Users should call `PetscOptionsLeftRestore()` to free the memory allocated in this routine
1978811af0c4SBarry Smith 
1979811af0c4SBarry Smith   The value returned may be different on each process and depends on which options have been processed
19801c9f3c13SBarry Smith   on the given process
19812d747510SLisandro Dalcin 
1982db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`
19832d747510SLisandro Dalcin @*/
1984d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftGet(PetscOptions options, PetscInt *N, char **names[], char **values[])
1985d71ae5a4SJacob Faibussowitsch {
19862d747510SLisandro Dalcin   PetscInt i, n;
19872d747510SLisandro Dalcin 
19882d747510SLisandro Dalcin   PetscFunctionBegin;
19894f572ea9SToby Isaac   if (N) PetscAssertPointer(N, 2);
19904f572ea9SToby Isaac   if (names) PetscAssertPointer(names, 3);
19914f572ea9SToby Isaac   if (values) PetscAssertPointer(values, 4);
19922d747510SLisandro Dalcin   options = options ? options : defaultoptions;
19932d747510SLisandro Dalcin 
19942d747510SLisandro Dalcin   /* The number of unused PETSc options */
19952d747510SLisandro Dalcin   n = 0;
19962d747510SLisandro Dalcin   for (i = 0; i < options->N; i++) {
1997660278c0SBarry Smith     if (PetscCIOption(options->names[i])) continue;
19982d747510SLisandro Dalcin     if (!options->used[i]) n++;
19992d747510SLisandro Dalcin   }
2000ad540459SPierre Jolivet   if (N) *N = n;
20019566063dSJacob Faibussowitsch   if (names) PetscCall(PetscMalloc1(n, names));
20029566063dSJacob Faibussowitsch   if (values) PetscCall(PetscMalloc1(n, values));
20032d747510SLisandro Dalcin 
20042d747510SLisandro Dalcin   n = 0;
20052d747510SLisandro Dalcin   if (names || values) {
20062d747510SLisandro Dalcin     for (i = 0; i < options->N; i++) {
20072d747510SLisandro Dalcin       if (!options->used[i]) {
2008660278c0SBarry Smith         if (PetscCIOption(options->names[i])) continue;
20092d747510SLisandro Dalcin         if (names) (*names)[n] = options->names[i];
20102d747510SLisandro Dalcin         if (values) (*values)[n] = options->values[i];
20112d747510SLisandro Dalcin         n++;
20122d747510SLisandro Dalcin       }
20132d747510SLisandro Dalcin     }
20142d747510SLisandro Dalcin   }
20153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20162d747510SLisandro Dalcin }
20172d747510SLisandro Dalcin 
20182d747510SLisandro Dalcin /*@C
2019811af0c4SBarry Smith   PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using `PetscOptionsLeftGet()`.
20202d747510SLisandro Dalcin 
20212d747510SLisandro Dalcin   Not Collective
20222d747510SLisandro Dalcin 
2023d8d19677SJose E. Roman   Input Parameters:
202420f4b53cSBarry Smith + options - options database, use `NULL` for default global database
202510450e9eSJacob Faibussowitsch . N       - count of options not used
20262d747510SLisandro Dalcin . names   - names of options not used
2027a2b725a8SWilliam Gropp - values  - values of options not used
20282d747510SLisandro Dalcin 
20292d747510SLisandro Dalcin   Level: advanced
20302d747510SLisandro Dalcin 
203110450e9eSJacob Faibussowitsch   Notes:
203210450e9eSJacob Faibussowitsch   The user should pass the same pointer to `N` as they did when calling `PetscOptionsLeftGet()`
203310450e9eSJacob Faibussowitsch 
2034db781477SPatrick Sanan .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`, `PetscOptionsLeftGet()`
20352d747510SLisandro Dalcin @*/
2036d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsLeftRestore(PetscOptions options, PetscInt *N, char **names[], char **values[])
2037d71ae5a4SJacob Faibussowitsch {
20382d747510SLisandro Dalcin   PetscFunctionBegin;
203910450e9eSJacob Faibussowitsch   (void)options;
20404f572ea9SToby Isaac   if (N) PetscAssertPointer(N, 2);
20414f572ea9SToby Isaac   if (names) PetscAssertPointer(names, 3);
20424f572ea9SToby Isaac   if (values) PetscAssertPointer(values, 4);
2043ad540459SPierre Jolivet   if (N) *N = 0;
20449566063dSJacob Faibussowitsch   if (names) PetscCall(PetscFree(*names));
20459566063dSJacob Faibussowitsch   if (values) PetscCall(PetscFree(*values));
20463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20472d747510SLisandro Dalcin }
20482d747510SLisandro Dalcin 
20492d747510SLisandro Dalcin /*@C
2050811af0c4SBarry Smith   PetscOptionsMonitorDefault - Print all options set value events using the supplied `PetscViewer`.
20512d747510SLisandro Dalcin 
2052c3339decSBarry Smith   Logically Collective
20532d747510SLisandro Dalcin 
20542d747510SLisandro Dalcin   Input Parameters:
20552d747510SLisandro Dalcin + name   - option name string
20562d747510SLisandro Dalcin . value  - option value string
20579355ec05SMatthew G. Knepley . source - The source for the option
205820f4b53cSBarry Smith - ctx    - a `PETSCVIEWERASCII` or `NULL`
20592d747510SLisandro Dalcin 
20602d747510SLisandro Dalcin   Level: intermediate
20612d747510SLisandro Dalcin 
20629666a313SBarry Smith   Notes:
206320f4b53cSBarry Smith   If ctx is `NULL`, `PetscPrintf()` is used.
20649314d9b7SBarry Smith   The first MPI process in the `PetscViewer` viewer actually prints the values, other
20651c9f3c13SBarry Smith   processes may have different values set
20661c9f3c13SBarry Smith 
2067811af0c4SBarry Smith   If `PetscCIEnabled` then do not print the test harness options
2068660278c0SBarry Smith 
2069db781477SPatrick Sanan .seealso: `PetscOptionsMonitorSet()`
20702d747510SLisandro Dalcin @*/
20719355ec05SMatthew G. Knepley PetscErrorCode PetscOptionsMonitorDefault(const char name[], const char value[], PetscOptionSource source, void *ctx)
2072d71ae5a4SJacob Faibussowitsch {
20732d747510SLisandro Dalcin   PetscFunctionBegin;
20743ba16761SJacob Faibussowitsch   if (PetscCIOption(name)) PetscFunctionReturn(PETSC_SUCCESS);
2075660278c0SBarry Smith 
20769060e2f9SVaclav Hapla   if (ctx) {
20779060e2f9SVaclav Hapla     PetscViewer viewer = (PetscViewer)ctx;
20782d747510SLisandro Dalcin     if (!value) {
20799566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Removing option: %s\n", name));
20802d747510SLisandro Dalcin     } else if (!value[0]) {
20819355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
20822d747510SLisandro Dalcin     } else {
20839355ec05SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
20842d747510SLisandro Dalcin     }
20859060e2f9SVaclav Hapla   } else {
20869060e2f9SVaclav Hapla     MPI_Comm comm = PETSC_COMM_WORLD;
20879060e2f9SVaclav Hapla     if (!value) {
20889566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Removing option: %s\n", name));
20899060e2f9SVaclav Hapla     } else if (!value[0]) {
20909355ec05SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
20919060e2f9SVaclav Hapla     } else {
2092aaa8cc7dSPierre Jolivet       PetscCall(PetscPrintf(comm, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
20939060e2f9SVaclav Hapla     }
20949060e2f9SVaclav Hapla   }
20953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20962d747510SLisandro Dalcin }
20972d747510SLisandro Dalcin 
20982d747510SLisandro Dalcin /*@C
20992d747510SLisandro Dalcin   PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
21002d747510SLisandro Dalcin   modified the PETSc options database.
21012d747510SLisandro Dalcin 
21022d747510SLisandro Dalcin   Not Collective
21032d747510SLisandro Dalcin 
21042d747510SLisandro Dalcin   Input Parameters:
210520f4b53cSBarry Smith + monitor        - pointer to function (if this is `NULL`, it turns off monitoring
210610450e9eSJacob Faibussowitsch . mctx           - [optional] context for private data for the monitor routine (use `NULL` if
210710450e9eSJacob Faibussowitsch                    no context is desired)
210810450e9eSJacob Faibussowitsch - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
21092d747510SLisandro Dalcin 
211010450e9eSJacob Faibussowitsch   Calling sequence of `monitor`:
21112d747510SLisandro Dalcin + name   - option name string
2112*432b765aSRené Chenard . value  - option value string, a value of `NULL` indicates the option is being removed from the database. A value
2113*432b765aSRené Chenard            of "" indicates the option is in the database but has no value.
21149355ec05SMatthew G. Knepley . source - option source
2115811af0c4SBarry Smith - mctx   - optional monitoring context, as set by `PetscOptionsMonitorSet()`
21162d747510SLisandro Dalcin 
211710450e9eSJacob Faibussowitsch   Calling sequence of `monitordestroy`:
211810450e9eSJacob Faibussowitsch . mctx - [optional] pointer to context to destroy with
21192d747510SLisandro Dalcin 
2120*432b765aSRené Chenard   Options Database Keys:
2121*432b765aSRené Chenard + -options_monitor <viewer> - turn on default monitoring
2122*432b765aSRené Chenard - -options_monitor_cancel   - turn off any option monitors except the default monitor obtained with `-options_monitor`
2123*432b765aSRené Chenard 
212420f4b53cSBarry Smith   Level: intermediate
212520f4b53cSBarry Smith 
21262d747510SLisandro Dalcin   Notes:
212710450e9eSJacob Faibussowitsch   See `PetscInitialize()` for options related to option database monitoring.
212810450e9eSJacob Faibussowitsch 
2129*432b765aSRené Chenard   The default is to do no monitoring.  To print the name and value of options
2130811af0c4SBarry Smith   being inserted into the database, use `PetscOptionsMonitorDefault()` as the monitoring routine,
2131*432b765aSRené Chenard   with a `NULL` monitoring context. Or use the option `-options_monitor` <viewer>.
21322d747510SLisandro Dalcin 
21332d747510SLisandro Dalcin   Several different monitoring routines may be set by calling
2134811af0c4SBarry Smith   `PetscOptionsMonitorSet()` multiple times; all will be called in the
21352d747510SLisandro Dalcin   order in which they were set.
21362d747510SLisandro Dalcin 
2137db781477SPatrick Sanan .seealso: `PetscOptionsMonitorDefault()`, `PetscInitialize()`
21382d747510SLisandro Dalcin @*/
213910450e9eSJacob Faibussowitsch PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], PetscOptionSource source, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx))
2140d71ae5a4SJacob Faibussowitsch {
21412d747510SLisandro Dalcin   PetscOptions options = defaultoptions;
21422d747510SLisandro Dalcin 
21432d747510SLisandro Dalcin   PetscFunctionBegin;
21443ba16761SJacob Faibussowitsch   if (options->monitorCancel) PetscFunctionReturn(PETSC_SUCCESS);
214508401ef6SPierre Jolivet   PetscCheck(options->numbermonitors < MAXOPTIONSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many PetscOptions monitors set");
21462d747510SLisandro Dalcin   options->monitor[options->numbermonitors]          = monitor;
21472d747510SLisandro Dalcin   options->monitordestroy[options->numbermonitors]   = monitordestroy;
21482d747510SLisandro Dalcin   options->monitorcontext[options->numbermonitors++] = (void *)mctx;
21493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21502d747510SLisandro Dalcin }
21512d747510SLisandro Dalcin 
21522d747510SLisandro Dalcin /*
21532d747510SLisandro Dalcin    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
215463fe8743SVaclav Hapla      Empty string is considered as true.
21552d747510SLisandro Dalcin */
2156d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToBool(const char value[], PetscBool *a)
2157d71ae5a4SJacob Faibussowitsch {
21582d747510SLisandro Dalcin   PetscBool istrue, isfalse;
21592d747510SLisandro Dalcin   size_t    len;
21602d747510SLisandro Dalcin 
21612d747510SLisandro Dalcin   PetscFunctionBegin;
216263fe8743SVaclav Hapla   /* PetscStrlen() returns 0 for NULL or "" */
21639566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(value, &len));
21649371c9d4SSatish Balay   if (!len) {
21659371c9d4SSatish Balay     *a = PETSC_TRUE;
21663ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21679371c9d4SSatish Balay   }
21689566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "TRUE", &istrue));
21699371c9d4SSatish Balay   if (istrue) {
21709371c9d4SSatish Balay     *a = PETSC_TRUE;
21713ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21729371c9d4SSatish Balay   }
21739566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "YES", &istrue));
21749371c9d4SSatish Balay   if (istrue) {
21759371c9d4SSatish Balay     *a = PETSC_TRUE;
21763ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21779371c9d4SSatish Balay   }
21789566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "1", &istrue));
21799371c9d4SSatish Balay   if (istrue) {
21809371c9d4SSatish Balay     *a = PETSC_TRUE;
21813ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21829371c9d4SSatish Balay   }
21839566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "on", &istrue));
21849371c9d4SSatish Balay   if (istrue) {
21859371c9d4SSatish Balay     *a = PETSC_TRUE;
21863ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21879371c9d4SSatish Balay   }
21889566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "FALSE", &isfalse));
21899371c9d4SSatish Balay   if (isfalse) {
21909371c9d4SSatish Balay     *a = PETSC_FALSE;
21913ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21929371c9d4SSatish Balay   }
21939566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "NO", &isfalse));
21949371c9d4SSatish Balay   if (isfalse) {
21959371c9d4SSatish Balay     *a = PETSC_FALSE;
21963ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
21979371c9d4SSatish Balay   }
21989566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "0", &isfalse));
21999371c9d4SSatish Balay   if (isfalse) {
22009371c9d4SSatish Balay     *a = PETSC_FALSE;
22013ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22029371c9d4SSatish Balay   }
22039566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(value, "off", &isfalse));
22049371c9d4SSatish Balay   if (isfalse) {
22059371c9d4SSatish Balay     *a = PETSC_FALSE;
22063ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
22079371c9d4SSatish Balay   }
220898921bdaSJacob Faibussowitsch   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown logical value: %s", value);
22092d747510SLisandro Dalcin }
22102d747510SLisandro Dalcin 
22112d747510SLisandro Dalcin /*
22122d747510SLisandro Dalcin    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
22132d747510SLisandro Dalcin */
2214d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToInt(const char name[], PetscInt *a)
2215d71ae5a4SJacob Faibussowitsch {
22162d747510SLisandro Dalcin   size_t    len;
22172d747510SLisandro Dalcin   PetscBool decide, tdefault, mouse;
22182d747510SLisandro Dalcin 
22192d747510SLisandro Dalcin   PetscFunctionBegin;
22209566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
22215f80ce2aSJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
22222d747510SLisandro Dalcin 
22239566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &tdefault));
222448a46eb9SPierre Jolivet   if (!tdefault) PetscCall(PetscStrcasecmp(name, "DEFAULT", &tdefault));
22259566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &decide));
222648a46eb9SPierre Jolivet   if (!decide) PetscCall(PetscStrcasecmp(name, "DECIDE", &decide));
22279566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "mouse", &mouse));
22282d747510SLisandro Dalcin 
22292d747510SLisandro Dalcin   if (tdefault) *a = PETSC_DEFAULT;
22302d747510SLisandro Dalcin   else if (decide) *a = PETSC_DECIDE;
22312d747510SLisandro Dalcin   else if (mouse) *a = -1;
22322d747510SLisandro Dalcin   else {
22332d747510SLisandro Dalcin     char *endptr;
22342d747510SLisandro Dalcin     long  strtolval;
22352d747510SLisandro Dalcin 
22362d747510SLisandro Dalcin     strtolval = strtol(name, &endptr, 10);
2237cc73adaaSBarry 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);
22382d747510SLisandro Dalcin 
22392d747510SLisandro Dalcin #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
22402d747510SLisandro Dalcin     (void)strtolval;
22412d747510SLisandro Dalcin     *a = atoll(name);
22422d747510SLisandro Dalcin #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
22432d747510SLisandro Dalcin     (void)strtolval;
22442d747510SLisandro Dalcin     *a = _atoi64(name);
22452d747510SLisandro Dalcin #else
22462d747510SLisandro Dalcin     *a = (PetscInt)strtolval;
22472d747510SLisandro Dalcin #endif
22482d747510SLisandro Dalcin   }
22493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22502d747510SLisandro Dalcin }
22512d747510SLisandro Dalcin 
22522d747510SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
22532d747510SLisandro Dalcin   #include <quadmath.h>
22542d747510SLisandro Dalcin #endif
22552d747510SLisandro Dalcin 
2256d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscStrtod(const char name[], PetscReal *a, char **endptr)
2257d71ae5a4SJacob Faibussowitsch {
22582d747510SLisandro Dalcin   PetscFunctionBegin;
22592d747510SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
22602d747510SLisandro Dalcin   *a = strtoflt128(name, endptr);
22612d747510SLisandro Dalcin #else
22622d747510SLisandro Dalcin   *a = (PetscReal)strtod(name, endptr);
22632d747510SLisandro Dalcin #endif
22643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22652d747510SLisandro Dalcin }
22662d747510SLisandro Dalcin 
2267d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscStrtoz(const char name[], PetscScalar *a, char **endptr, PetscBool *isImaginary)
2268d71ae5a4SJacob Faibussowitsch {
22692d747510SLisandro Dalcin   PetscBool hasi = PETSC_FALSE;
22702d747510SLisandro Dalcin   char     *ptr;
22712d747510SLisandro Dalcin   PetscReal strtoval;
22722d747510SLisandro Dalcin 
22732d747510SLisandro Dalcin   PetscFunctionBegin;
22749566063dSJacob Faibussowitsch   PetscCall(PetscStrtod(name, &strtoval, &ptr));
22752d747510SLisandro Dalcin   if (ptr == name) {
22762d747510SLisandro Dalcin     strtoval = 1.;
22772d747510SLisandro Dalcin     hasi     = PETSC_TRUE;
22782d747510SLisandro Dalcin     if (name[0] == 'i') {
22792d747510SLisandro Dalcin       ptr++;
22802d747510SLisandro Dalcin     } else if (name[0] == '+' && name[1] == 'i') {
22812d747510SLisandro Dalcin       ptr += 2;
22822d747510SLisandro Dalcin     } else if (name[0] == '-' && name[1] == 'i') {
22832d747510SLisandro Dalcin       strtoval = -1.;
22842d747510SLisandro Dalcin       ptr += 2;
22852d747510SLisandro Dalcin     }
22862d747510SLisandro Dalcin   } else if (*ptr == 'i') {
22872d747510SLisandro Dalcin     hasi = PETSC_TRUE;
22882d747510SLisandro Dalcin     ptr++;
22892d747510SLisandro Dalcin   }
22902d747510SLisandro Dalcin   *endptr      = ptr;
22912d747510SLisandro Dalcin   *isImaginary = hasi;
22922d747510SLisandro Dalcin   if (hasi) {
22932d747510SLisandro Dalcin #if !defined(PETSC_USE_COMPLEX)
229498921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s contains imaginary but complex not supported ", name);
22952d747510SLisandro Dalcin #else
22962d747510SLisandro Dalcin     *a = PetscCMPLX(0., strtoval);
22972d747510SLisandro Dalcin #endif
22982d747510SLisandro Dalcin   } else {
22992d747510SLisandro Dalcin     *a = strtoval;
23002d747510SLisandro Dalcin   }
23013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23022d747510SLisandro Dalcin }
23032d747510SLisandro Dalcin 
23042d747510SLisandro Dalcin /*
23052d747510SLisandro Dalcin    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
23062d747510SLisandro Dalcin */
2307d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToReal(const char name[], PetscReal *a)
2308d71ae5a4SJacob Faibussowitsch {
23092d747510SLisandro Dalcin   size_t    len;
23102d747510SLisandro Dalcin   PetscBool match;
23112d747510SLisandro Dalcin   char     *endptr;
23122d747510SLisandro Dalcin 
23132d747510SLisandro Dalcin   PetscFunctionBegin;
23149566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
231528b400f6SJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "String of length zero has no numerical value");
23162d747510SLisandro Dalcin 
23179566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &match));
23189566063dSJacob Faibussowitsch   if (!match) PetscCall(PetscStrcasecmp(name, "DEFAULT", &match));
23199371c9d4SSatish Balay   if (match) {
23209371c9d4SSatish Balay     *a = PETSC_DEFAULT;
23213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
23229371c9d4SSatish Balay   }
23232d747510SLisandro Dalcin 
23249566063dSJacob Faibussowitsch   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &match));
23259566063dSJacob Faibussowitsch   if (!match) PetscCall(PetscStrcasecmp(name, "DECIDE", &match));
23269371c9d4SSatish Balay   if (match) {
23279371c9d4SSatish Balay     *a = PETSC_DECIDE;
23283ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
23299371c9d4SSatish Balay   }
23302d747510SLisandro Dalcin 
23319566063dSJacob Faibussowitsch   PetscCall(PetscStrtod(name, a, &endptr));
233239a651e2SJacob Faibussowitsch   PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value", name);
23333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23342d747510SLisandro Dalcin }
23352d747510SLisandro Dalcin 
2336d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringToScalar(const char name[], PetscScalar *a)
2337d71ae5a4SJacob Faibussowitsch {
23382d747510SLisandro Dalcin   PetscBool   imag1;
23392d747510SLisandro Dalcin   size_t      len;
23402d747510SLisandro Dalcin   PetscScalar val = 0.;
23412d747510SLisandro Dalcin   char       *ptr = NULL;
23422d747510SLisandro Dalcin 
23432d747510SLisandro Dalcin   PetscFunctionBegin;
23449566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
234528b400f6SJacob Faibussowitsch   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
23469566063dSJacob Faibussowitsch   PetscCall(PetscStrtoz(name, &val, &ptr, &imag1));
23472d747510SLisandro Dalcin #if defined(PETSC_USE_COMPLEX)
23482d747510SLisandro Dalcin   if ((size_t)(ptr - name) < len) {
23492d747510SLisandro Dalcin     PetscBool   imag2;
23502d747510SLisandro Dalcin     PetscScalar val2;
23512d747510SLisandro Dalcin 
23529566063dSJacob Faibussowitsch     PetscCall(PetscStrtoz(ptr, &val2, &ptr, &imag2));
235339a651e2SJacob Faibussowitsch     if (imag1) PetscCheck(imag2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s: must specify imaginary component second", name);
23542d747510SLisandro Dalcin     val = PetscCMPLX(PetscRealPart(val), PetscImaginaryPart(val2));
23552d747510SLisandro Dalcin   }
23562d747510SLisandro Dalcin #endif
235739a651e2SJacob Faibussowitsch   PetscCheck((size_t)(ptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value ", name);
23582d747510SLisandro Dalcin   *a = val;
23593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23602d747510SLisandro Dalcin }
23612d747510SLisandro Dalcin 
23622d747510SLisandro Dalcin /*@C
23632d747510SLisandro Dalcin   PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
23642d747510SLisandro Dalcin   option in the database.
2365e5c89e4eSSatish Balay 
2366e5c89e4eSSatish Balay   Not Collective
2367e5c89e4eSSatish Balay 
2368e5c89e4eSSatish Balay   Input Parameters:
236920f4b53cSBarry Smith + options - options database, use `NULL` for default global database
237020f4b53cSBarry Smith . pre     - the string to prepend to the name or `NULL`
2371e5c89e4eSSatish Balay - name    - the option one is seeking
2372e5c89e4eSSatish Balay 
2373d8d19677SJose E. Roman   Output Parameters:
23742d747510SLisandro Dalcin + ivalue - the logical value to return
2375811af0c4SBarry Smith - set    - `PETSC_TRUE`  if found, else `PETSC_FALSE`
2376e5c89e4eSSatish Balay 
2377e5c89e4eSSatish Balay   Level: beginner
2378e5c89e4eSSatish Balay 
237995452b02SPatrick Sanan   Notes:
2380811af0c4SBarry Smith   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
2381811af0c4SBarry Smith   FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
23822d747510SLisandro Dalcin 
23839314d9b7SBarry 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`
23849314d9b7SBarry Smith   is equivalent to `-requested_bool true`
23852d747510SLisandro Dalcin 
23869314d9b7SBarry Smith   If the user does not supply the option at all `ivalue` is NOT changed. Thus
23879314d9b7SBarry Smith   you should ALWAYS initialize `ivalue` if you access it without first checking that the `set` flag is true.
23882efd9cb1SBarry Smith 
2389db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2390db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2391db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2392c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2393db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2394db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2395e5c89e4eSSatish Balay @*/
2396d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetBool(PetscOptions options, const char pre[], const char name[], PetscBool *ivalue, PetscBool *set)
2397d71ae5a4SJacob Faibussowitsch {
23982d747510SLisandro Dalcin   const char *value;
2399ace3abfcSBarry Smith   PetscBool   flag;
2400e5c89e4eSSatish Balay 
2401e5c89e4eSSatish Balay   PetscFunctionBegin;
24024f572ea9SToby Isaac   PetscAssertPointer(name, 3);
24034f572ea9SToby Isaac   if (ivalue) PetscAssertPointer(ivalue, 4);
24049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2405e5c89e4eSSatish Balay   if (flag) {
240696ef3cdfSSatish Balay     if (set) *set = PETSC_TRUE;
24079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(value, &flag));
24082d747510SLisandro Dalcin     if (ivalue) *ivalue = flag;
2409e5c89e4eSSatish Balay   } else {
241096ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2411e5c89e4eSSatish Balay   }
24123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2413e5c89e4eSSatish Balay }
2414e5c89e4eSSatish Balay 
2415e5c89e4eSSatish Balay /*@C
2416e5c89e4eSSatish Balay   PetscOptionsGetEList - Puts a list of option values that a single one may be selected from
2417e5c89e4eSSatish Balay 
2418e5c89e4eSSatish Balay   Not Collective
2419e5c89e4eSSatish Balay 
2420e5c89e4eSSatish Balay   Input Parameters:
242120f4b53cSBarry Smith + options - options database, use `NULL` for default global database
242220f4b53cSBarry Smith . pre     - the string to prepend to the name or `NULL`
2423e5c89e4eSSatish Balay . opt     - option name
2424a264d7a6SBarry Smith . list    - the possible choices (one of these must be selected, anything else is invalid)
2425a2b725a8SWilliam Gropp - ntext   - number of choices
2426e5c89e4eSSatish Balay 
2427d8d19677SJose E. Roman   Output Parameters:
24282efd9cb1SBarry Smith + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2429811af0c4SBarry Smith - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
2430e5c89e4eSSatish Balay 
2431e5c89e4eSSatish Balay   Level: intermediate
2432e5c89e4eSSatish Balay 
243395452b02SPatrick Sanan   Notes:
24349314d9b7SBarry Smith   If the user does not supply the option `value` is NOT changed. Thus
24359314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is true.
24362efd9cb1SBarry Smith 
2437811af0c4SBarry Smith   See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList`
2438e5c89e4eSSatish Balay 
2439db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2440db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2441db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2442c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2443db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2444db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2445e5c89e4eSSatish Balay @*/
2446d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEList(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscInt ntext, PetscInt *value, PetscBool *set)
2447d71ae5a4SJacob Faibussowitsch {
244858b0ac4eSStefano Zampini   size_t    alen, len = 0, tlen = 0;
2449e5c89e4eSSatish Balay   char     *svalue;
2450ace3abfcSBarry Smith   PetscBool aset, flg = PETSC_FALSE;
2451e5c89e4eSSatish Balay   PetscInt  i;
2452e5c89e4eSSatish Balay 
2453e5c89e4eSSatish Balay   PetscFunctionBegin;
24544f572ea9SToby Isaac   PetscAssertPointer(opt, 3);
2455e5c89e4eSSatish Balay   for (i = 0; i < ntext; i++) {
24569566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(list[i], &alen));
2457e5c89e4eSSatish Balay     if (alen > len) len = alen;
245858b0ac4eSStefano Zampini     tlen += len + 1;
2459e5c89e4eSSatish Balay   }
2460e5c89e4eSSatish Balay   len += 5; /* a little extra space for user mistypes */
24619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &svalue));
24629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(options, pre, opt, svalue, len, &aset));
2463e5c89e4eSSatish Balay   if (aset) {
24649566063dSJacob Faibussowitsch     PetscCall(PetscEListFind(ntext, list, svalue, value, &flg));
246558b0ac4eSStefano Zampini     if (!flg) {
2466c6a7a370SJeremy L Thompson       char *avail;
246758b0ac4eSStefano Zampini 
24689566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(tlen, &avail));
2469c6a7a370SJeremy L Thompson       avail[0] = '\0';
247058b0ac4eSStefano Zampini       for (i = 0; i < ntext; i++) {
2471c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(avail, list[i], tlen));
2472c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(avail, " ", tlen));
247358b0ac4eSStefano Zampini       }
24749566063dSJacob Faibussowitsch       PetscCall(PetscStrtolower(avail));
247598921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown option %s for -%s%s. Available options: %s", svalue, pre ? pre : "", opt + 1, avail);
247658b0ac4eSStefano Zampini     }
2477fbedd5e0SJed Brown     if (set) *set = PETSC_TRUE;
2478a297a907SKarl Rupp   } else if (set) *set = PETSC_FALSE;
24799566063dSJacob Faibussowitsch   PetscCall(PetscFree(svalue));
24803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2481e5c89e4eSSatish Balay }
2482e5c89e4eSSatish Balay 
2483e5c89e4eSSatish Balay /*@C
2484e5c89e4eSSatish Balay   PetscOptionsGetEnum - Gets the enum value for a particular option in the database.
2485e5c89e4eSSatish Balay 
2486e5c89e4eSSatish Balay   Not Collective
2487e5c89e4eSSatish Balay 
2488e5c89e4eSSatish Balay   Input Parameters:
248920f4b53cSBarry Smith + options - options database, use `NULL` for default global database
249020f4b53cSBarry Smith . pre     - option prefix or `NULL`
2491e5c89e4eSSatish Balay . opt     - option name
24926b867d5aSJose E. Roman - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2493e5c89e4eSSatish Balay 
2494d8d19677SJose E. Roman   Output Parameters:
2495e5c89e4eSSatish Balay + value - the value to return
2496811af0c4SBarry Smith - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`
2497e5c89e4eSSatish Balay 
2498e5c89e4eSSatish Balay   Level: beginner
2499e5c89e4eSSatish Balay 
250095452b02SPatrick Sanan   Notes:
25019314d9b7SBarry Smith   If the user does not supply the option `value` is NOT changed. Thus
25029314d9b7SBarry Smith   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is true.
2503e5c89e4eSSatish Balay 
25049314d9b7SBarry Smith   `list` is usually something like `PCASMTypes` or some other predefined list of enum names
2505e5c89e4eSSatish Balay 
2506db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2507db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2508aec76313SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2509db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2510c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2511db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2512db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2513e5c89e4eSSatish Balay @*/
2514d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEnum(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscEnum *value, PetscBool *set)
2515d71ae5a4SJacob Faibussowitsch {
251669a24498SJed Brown   PetscInt  ntext = 0, tval;
2517ace3abfcSBarry Smith   PetscBool fset;
2518e5c89e4eSSatish Balay 
2519e5c89e4eSSatish Balay   PetscFunctionBegin;
25204f572ea9SToby Isaac   PetscAssertPointer(opt, 3);
2521ad540459SPierre 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");
252208401ef6SPierre Jolivet   PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix");
2523e5c89e4eSSatish Balay   ntext -= 3;
25249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEList(options, pre, opt, list, ntext, &tval, &fset));
252569a24498SJed Brown   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2526809ceb46SBarry Smith   if (fset) *value = (PetscEnum)tval;
2527809ceb46SBarry Smith   if (set) *set = fset;
25283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2529e5c89e4eSSatish Balay }
2530e5c89e4eSSatish Balay 
2531e5c89e4eSSatish Balay /*@C
25322d747510SLisandro Dalcin   PetscOptionsGetInt - Gets the integer value for a particular option in the database.
2533e5c89e4eSSatish Balay 
2534e5c89e4eSSatish Balay   Not Collective
2535e5c89e4eSSatish Balay 
2536e5c89e4eSSatish Balay   Input Parameters:
253720f4b53cSBarry Smith + options - options database, use `NULL` for default global database
253820f4b53cSBarry Smith . pre     - the string to prepend to the name or `NULL`
2539e5c89e4eSSatish Balay - name    - the option one is seeking
2540e5c89e4eSSatish Balay 
2541d8d19677SJose E. Roman   Output Parameters:
25422d747510SLisandro Dalcin + ivalue - the integer value to return
2543811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
2544e5c89e4eSSatish Balay 
2545e5c89e4eSSatish Balay   Level: beginner
2546e5c89e4eSSatish Balay 
2547e5c89e4eSSatish Balay   Notes:
25489314d9b7SBarry Smith   If the user does not supply the option `ivalue` is NOT changed. Thus
25499314d9b7SBarry Smith   you should ALWAYS initialize the `ivalue` if you access it without first checking that the `set` flag is true.
25505c07ccb8SBarry Smith 
2551db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2552db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2553aec76313SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2554db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2555c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2556db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2557db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2558e5c89e4eSSatish Balay @*/
2559d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetInt(PetscOptions options, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
2560d71ae5a4SJacob Faibussowitsch {
25612d747510SLisandro Dalcin   const char *value;
25622d747510SLisandro Dalcin   PetscBool   flag;
2563e5c89e4eSSatish Balay 
2564e5c89e4eSSatish Balay   PetscFunctionBegin;
25654f572ea9SToby Isaac   PetscAssertPointer(name, 3);
25664f572ea9SToby Isaac   PetscAssertPointer(ivalue, 4);
25679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2568e5c89e4eSSatish Balay   if (flag) {
256934a9cc2cSBarry Smith     if (!value) {
25702d747510SLisandro Dalcin       if (set) *set = PETSC_FALSE;
257134a9cc2cSBarry Smith     } else {
25722d747510SLisandro Dalcin       if (set) *set = PETSC_TRUE;
25739566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToInt(value, ivalue));
2574e5c89e4eSSatish Balay     }
2575e5c89e4eSSatish Balay   } else {
257696ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2577e5c89e4eSSatish Balay   }
25783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2579e5c89e4eSSatish Balay }
2580e5c89e4eSSatish Balay 
2581e2446a98SMatthew Knepley /*@C
2582e5c89e4eSSatish Balay   PetscOptionsGetReal - Gets the double precision value for a particular
2583e5c89e4eSSatish Balay   option in the database.
2584e5c89e4eSSatish Balay 
2585e5c89e4eSSatish Balay   Not Collective
2586e5c89e4eSSatish Balay 
2587e5c89e4eSSatish Balay   Input Parameters:
258820f4b53cSBarry Smith + options - options database, use `NULL` for default global database
258920f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
2590e5c89e4eSSatish Balay - name    - the option one is seeking
2591e5c89e4eSSatish Balay 
2592d8d19677SJose E. Roman   Output Parameters:
2593e5c89e4eSSatish Balay + dvalue - the double value to return
2594811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, `PETSC_FALSE` if not found
2595e5c89e4eSSatish Balay 
259620f4b53cSBarry Smith   Level: beginner
259720f4b53cSBarry Smith 
2598811af0c4SBarry Smith   Note:
25999314d9b7SBarry Smith   If the user does not supply the option `dvalue` is NOT changed. Thus
26009314d9b7SBarry Smith   you should ALWAYS initialize `dvalue` if you access it without first checking that the `set` flag is true.
2601e4974155SBarry Smith 
2602db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2603c2e3fba1SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2604db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2605c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2606db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2607db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2608e5c89e4eSSatish Balay @*/
2609d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetReal(PetscOptions options, const char pre[], const char name[], PetscReal *dvalue, PetscBool *set)
2610d71ae5a4SJacob Faibussowitsch {
26112d747510SLisandro Dalcin   const char *value;
2612ace3abfcSBarry Smith   PetscBool   flag;
2613e5c89e4eSSatish Balay 
2614e5c89e4eSSatish Balay   PetscFunctionBegin;
26154f572ea9SToby Isaac   PetscAssertPointer(name, 3);
26164f572ea9SToby Isaac   PetscAssertPointer(dvalue, 4);
26179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2618e5c89e4eSSatish Balay   if (flag) {
2619a297a907SKarl Rupp     if (!value) {
2620a297a907SKarl Rupp       if (set) *set = PETSC_FALSE;
2621a297a907SKarl Rupp     } else {
2622a297a907SKarl Rupp       if (set) *set = PETSC_TRUE;
26239566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToReal(value, dvalue));
2624a297a907SKarl Rupp     }
2625e5c89e4eSSatish Balay   } else {
262696ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2627e5c89e4eSSatish Balay   }
26283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2629e5c89e4eSSatish Balay }
2630e5c89e4eSSatish Balay 
2631e5c89e4eSSatish Balay /*@C
2632e5c89e4eSSatish Balay   PetscOptionsGetScalar - Gets the scalar value for a particular
2633e5c89e4eSSatish Balay   option in the database.
2634e5c89e4eSSatish Balay 
2635e5c89e4eSSatish Balay   Not Collective
2636e5c89e4eSSatish Balay 
2637e5c89e4eSSatish Balay   Input Parameters:
263820f4b53cSBarry Smith + options - options database, use `NULL` for default global database
263920f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
2640e5c89e4eSSatish Balay - name    - the option one is seeking
2641e5c89e4eSSatish Balay 
2642d8d19677SJose E. Roman   Output Parameters:
26439314d9b7SBarry Smith + dvalue - the scalar value to return
2644811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
2645e5c89e4eSSatish Balay 
2646e5c89e4eSSatish Balay   Level: beginner
2647e5c89e4eSSatish Balay 
264810450e9eSJacob Faibussowitsch   Example Usage:
2649eb4ae41dSBarry Smith   A complex number 2+3i must be specified with NO spaces
2650e5c89e4eSSatish Balay 
2651811af0c4SBarry Smith   Note:
26529314d9b7SBarry Smith   If the user does not supply the option `dvalue` is NOT changed. Thus
26539314d9b7SBarry Smith   you should ALWAYS initialize `dvalue` if you access it without first checking if the `set` flag is true.
2654e4974155SBarry Smith 
2655db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2656db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2657db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2658c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2659db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2660db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2661e5c89e4eSSatish Balay @*/
2662d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetScalar(PetscOptions options, const char pre[], const char name[], PetscScalar *dvalue, PetscBool *set)
2663d71ae5a4SJacob Faibussowitsch {
26642d747510SLisandro Dalcin   const char *value;
2665ace3abfcSBarry Smith   PetscBool   flag;
2666e5c89e4eSSatish Balay 
2667e5c89e4eSSatish Balay   PetscFunctionBegin;
26684f572ea9SToby Isaac   PetscAssertPointer(name, 3);
26694f572ea9SToby Isaac   PetscAssertPointer(dvalue, 4);
26709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2671e5c89e4eSSatish Balay   if (flag) {
2672e5c89e4eSSatish Balay     if (!value) {
267396ef3cdfSSatish Balay       if (set) *set = PETSC_FALSE;
2674e5c89e4eSSatish Balay     } else {
2675e5c89e4eSSatish Balay #if !defined(PETSC_USE_COMPLEX)
26769566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToReal(value, dvalue));
2677e5c89e4eSSatish Balay #else
26789566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToScalar(value, dvalue));
2679e5c89e4eSSatish Balay #endif
268096ef3cdfSSatish Balay       if (set) *set = PETSC_TRUE;
2681e5c89e4eSSatish Balay     }
2682e5c89e4eSSatish Balay   } else { /* flag */
268396ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2684e5c89e4eSSatish Balay   }
26853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2686e5c89e4eSSatish Balay }
2687e5c89e4eSSatish Balay 
2688e5c89e4eSSatish Balay /*@C
2689e5c89e4eSSatish Balay   PetscOptionsGetString - Gets the string value for a particular option in
2690e5c89e4eSSatish Balay   the database.
2691e5c89e4eSSatish Balay 
2692e5c89e4eSSatish Balay   Not Collective
2693e5c89e4eSSatish Balay 
2694e5c89e4eSSatish Balay   Input Parameters:
269520f4b53cSBarry Smith + options - options database, use `NULL` for default global database
269620f4b53cSBarry Smith . pre     - string to prepend to name or `NULL`
2697e5c89e4eSSatish Balay . name    - the option one is seeking
2698bcbf2dc5SJed Brown - len     - maximum length of the string including null termination
2699e5c89e4eSSatish Balay 
2700e5c89e4eSSatish Balay   Output Parameters:
2701e5c89e4eSSatish Balay + string - location to copy string
2702811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
2703e5c89e4eSSatish Balay 
2704e5c89e4eSSatish Balay   Level: beginner
2705e5c89e4eSSatish Balay 
270620f4b53cSBarry Smith   Note:
27079314d9b7SBarry 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`
270820f4b53cSBarry Smith 
27099314d9b7SBarry Smith   If the user does not use the option then `string` is not changed. Thus
27109314d9b7SBarry Smith   you should ALWAYS initialize `string` if you access it without first checking that the `set` flag is true.
271120f4b53cSBarry Smith 
2712aec76313SJacob Faibussowitsch   Fortran Notes:
2713e5c89e4eSSatish Balay   The Fortran interface is slightly different from the C/C++
2714e5c89e4eSSatish Balay   interface (len is not used).  Sample usage in Fortran follows
2715e5c89e4eSSatish Balay .vb
2716e5c89e4eSSatish Balay       character *20    string
271793e6ba5cSBarry Smith       PetscErrorCode   ierr
271893e6ba5cSBarry Smith       PetscBool        set
27191b266c99SBarry Smith       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2720e5c89e4eSSatish Balay .ve
2721e5c89e4eSSatish Balay 
2722db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2723db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2724db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2725c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2726db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2727db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
2728e5c89e4eSSatish Balay @*/
2729d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetString(PetscOptions options, const char pre[], const char name[], char string[], size_t len, PetscBool *set)
2730d71ae5a4SJacob Faibussowitsch {
27312d747510SLisandro Dalcin   const char *value;
2732ace3abfcSBarry Smith   PetscBool   flag;
2733e5c89e4eSSatish Balay 
2734e5c89e4eSSatish Balay   PetscFunctionBegin;
27354f572ea9SToby Isaac   PetscAssertPointer(name, 3);
27364f572ea9SToby Isaac   PetscAssertPointer(string, 4);
27379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2738e5c89e4eSSatish Balay   if (!flag) {
273996ef3cdfSSatish Balay     if (set) *set = PETSC_FALSE;
2740e5c89e4eSSatish Balay   } else {
274196ef3cdfSSatish Balay     if (set) *set = PETSC_TRUE;
27429566063dSJacob Faibussowitsch     if (value) PetscCall(PetscStrncpy(string, value, len));
27439566063dSJacob Faibussowitsch     else PetscCall(PetscArrayzero(string, len));
2744e5c89e4eSSatish Balay   }
27453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2746e5c89e4eSSatish Balay }
2747e5c89e4eSSatish Balay 
27482d747510SLisandro Dalcin /*@C
27492d747510SLisandro Dalcin   PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2750f1a722f8SMatthew G. Knepley   option in the database.  The values must be separated with commas with no intervening spaces.
27512d747510SLisandro Dalcin 
27522d747510SLisandro Dalcin   Not Collective
27532d747510SLisandro Dalcin 
27542d747510SLisandro Dalcin   Input Parameters:
275520f4b53cSBarry Smith + options - options database, use `NULL` for default global database
275620f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
27576b867d5aSJose E. Roman - name    - the option one is seeking
27586b867d5aSJose E. Roman 
2759d8d19677SJose E. Roman   Output Parameters:
27609314d9b7SBarry Smith + dvalue - the Boolean values to return
2761f1a722f8SMatthew G. Knepley . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2762811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
27632d747510SLisandro Dalcin 
27642d747510SLisandro Dalcin   Level: beginner
27652d747510SLisandro Dalcin 
2766811af0c4SBarry Smith   Note:
276720f4b53cSBarry Smith   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`. FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
27682d747510SLisandro Dalcin 
2769db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2770db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2771db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2772c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2773db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2774db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
27752d747510SLisandro Dalcin @*/
2776d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options, const char pre[], const char name[], PetscBool dvalue[], PetscInt *nmax, PetscBool *set)
2777d71ae5a4SJacob Faibussowitsch {
27782d747510SLisandro Dalcin   const char *svalue;
27792d747510SLisandro Dalcin   char       *value;
27802d747510SLisandro Dalcin   PetscInt    n = 0;
27812d747510SLisandro Dalcin   PetscBool   flag;
27822d747510SLisandro Dalcin   PetscToken  token;
27832d747510SLisandro Dalcin 
27842d747510SLisandro Dalcin   PetscFunctionBegin;
27854f572ea9SToby Isaac   PetscAssertPointer(name, 3);
27864f572ea9SToby Isaac   PetscAssertPointer(dvalue, 4);
27874f572ea9SToby Isaac   PetscAssertPointer(nmax, 5);
27882d747510SLisandro Dalcin 
27899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
27909371c9d4SSatish Balay   if (!flag || !svalue) {
27919371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
27929371c9d4SSatish Balay     *nmax = 0;
27933ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
27949371c9d4SSatish Balay   }
27952d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
27969566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
27979566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
27982d747510SLisandro Dalcin   while (value && n < *nmax) {
27999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToBool(value, dvalue));
28009566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
28012d747510SLisandro Dalcin     dvalue++;
28022d747510SLisandro Dalcin     n++;
28032d747510SLisandro Dalcin   }
28049566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
28052d747510SLisandro Dalcin   *nmax = n;
28063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28072d747510SLisandro Dalcin }
28082d747510SLisandro Dalcin 
28092d747510SLisandro Dalcin /*@C
28102d747510SLisandro Dalcin   PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.
28112d747510SLisandro Dalcin 
28122d747510SLisandro Dalcin   Not Collective
28132d747510SLisandro Dalcin 
28142d747510SLisandro Dalcin   Input Parameters:
281520f4b53cSBarry Smith + options - options database, use `NULL` for default global database
281620f4b53cSBarry Smith . pre     - option prefix or `NULL`
28172d747510SLisandro Dalcin . name    - option name
28186b867d5aSJose E. Roman - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
28196b867d5aSJose E. Roman 
28202d747510SLisandro Dalcin   Output Parameters:
28212d747510SLisandro Dalcin + ivalue - the  enum values to return
2822f1a722f8SMatthew G. Knepley . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2823811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
28242d747510SLisandro Dalcin 
28252d747510SLisandro Dalcin   Level: beginner
28262d747510SLisandro Dalcin 
28272d747510SLisandro Dalcin   Notes:
28289314d9b7SBarry Smith   The array must be passed as a comma separated list with no spaces between the items.
28292d747510SLisandro Dalcin 
28309314d9b7SBarry Smith   `list` is usually something like `PCASMTypes` or some other predefined list of enum names.
28312d747510SLisandro Dalcin 
2832db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2833db781477SPatrick Sanan           `PetscOptionsGetEnum()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2834aec76313SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsName()`,
2835c2e3fba1SPatrick Sanan           `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, `PetscOptionsStringArray()`, `PetscOptionsRealArray()`,
2836db781477SPatrick Sanan           `PetscOptionsScalar()`, `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2837db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
28382d747510SLisandro Dalcin @*/
2839d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options, const char pre[], const char name[], const char *const *list, PetscEnum ivalue[], PetscInt *nmax, PetscBool *set)
2840d71ae5a4SJacob Faibussowitsch {
28412d747510SLisandro Dalcin   const char *svalue;
28422d747510SLisandro Dalcin   char       *value;
28432d747510SLisandro Dalcin   PetscInt    n = 0;
28442d747510SLisandro Dalcin   PetscEnum   evalue;
28452d747510SLisandro Dalcin   PetscBool   flag;
28462d747510SLisandro Dalcin   PetscToken  token;
28472d747510SLisandro Dalcin 
28482d747510SLisandro Dalcin   PetscFunctionBegin;
28494f572ea9SToby Isaac   PetscAssertPointer(name, 3);
28504f572ea9SToby Isaac   PetscAssertPointer(list, 4);
28514f572ea9SToby Isaac   PetscAssertPointer(ivalue, 5);
28524f572ea9SToby Isaac   PetscAssertPointer(nmax, 6);
28532d747510SLisandro Dalcin 
28549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
28559371c9d4SSatish Balay   if (!flag || !svalue) {
28569371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
28579371c9d4SSatish Balay     *nmax = 0;
28583ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
28599371c9d4SSatish Balay   }
28602d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
28619566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
28629566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
28632d747510SLisandro Dalcin   while (value && n < *nmax) {
28649566063dSJacob Faibussowitsch     PetscCall(PetscEnumFind(list, value, &evalue, &flag));
286528b400f6SJacob Faibussowitsch     PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown enum value '%s' for -%s%s", svalue, pre ? pre : "", name + 1);
28662d747510SLisandro Dalcin     ivalue[n++] = evalue;
28679566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
28682d747510SLisandro Dalcin   }
28699566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
28702d747510SLisandro Dalcin   *nmax = n;
28713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28722d747510SLisandro Dalcin }
28732d747510SLisandro Dalcin 
28742d747510SLisandro Dalcin /*@C
2875f1a722f8SMatthew G. Knepley   PetscOptionsGetIntArray - Gets an array of integer values for a particular option in the database.
28762d747510SLisandro Dalcin 
28772d747510SLisandro Dalcin   Not Collective
28782d747510SLisandro Dalcin 
28792d747510SLisandro Dalcin   Input Parameters:
288020f4b53cSBarry Smith + options - options database, use `NULL` for default global database
288120f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
28826b867d5aSJose E. Roman - name    - the option one is seeking
28836b867d5aSJose E. Roman 
2884d8d19677SJose E. Roman   Output Parameters:
28852d747510SLisandro Dalcin + ivalue - the integer values to return
2886f1a722f8SMatthew G. Knepley . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2887811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
28882d747510SLisandro Dalcin 
28892d747510SLisandro Dalcin   Level: beginner
28902d747510SLisandro Dalcin 
28912d747510SLisandro Dalcin   Notes:
28922d747510SLisandro Dalcin   The array can be passed as
2893811af0c4SBarry Smith +  a comma separated list -                                 0,1,2,3,4,5,6,7
2894811af0c4SBarry Smith .  a range (start\-end+1) -                                 0-8
2895811af0c4SBarry Smith .  a range with given increment (start\-end+1:inc) -        0-7:2
2896811af0c4SBarry Smith -  a combination of values and ranges separated by commas - 0,1-8,8-15:2
28972d747510SLisandro Dalcin 
28982d747510SLisandro Dalcin   There must be no intervening spaces between the values.
28992d747510SLisandro Dalcin 
2900db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2901db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2902db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2903c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2904db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2905db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
29062d747510SLisandro Dalcin @*/
2907d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetIntArray(PetscOptions options, const char pre[], const char name[], PetscInt ivalue[], PetscInt *nmax, PetscBool *set)
2908d71ae5a4SJacob Faibussowitsch {
29092d747510SLisandro Dalcin   const char *svalue;
29102d747510SLisandro Dalcin   char       *value;
29112d747510SLisandro Dalcin   PetscInt    n = 0, i, j, start, end, inc, nvalues;
29122d747510SLisandro Dalcin   size_t      len;
29132d747510SLisandro Dalcin   PetscBool   flag, foundrange;
29142d747510SLisandro Dalcin   PetscToken  token;
29152d747510SLisandro Dalcin 
29162d747510SLisandro Dalcin   PetscFunctionBegin;
29174f572ea9SToby Isaac   PetscAssertPointer(name, 3);
29184f572ea9SToby Isaac   PetscAssertPointer(ivalue, 4);
29194f572ea9SToby Isaac   PetscAssertPointer(nmax, 5);
29202d747510SLisandro Dalcin 
29219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
29229371c9d4SSatish Balay   if (!flag || !svalue) {
29239371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
29249371c9d4SSatish Balay     *nmax = 0;
29253ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
29269371c9d4SSatish Balay   }
29272d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
29289566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
29299566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
29302d747510SLisandro Dalcin   while (value && n < *nmax) {
29312d747510SLisandro Dalcin     /* look for form  d-D where d and D are integers */
29322d747510SLisandro Dalcin     foundrange = PETSC_FALSE;
29339566063dSJacob Faibussowitsch     PetscCall(PetscStrlen(value, &len));
29342d747510SLisandro Dalcin     if (value[0] == '-') i = 2;
29352d747510SLisandro Dalcin     else i = 1;
29362d747510SLisandro Dalcin     for (; i < (int)len; i++) {
29372d747510SLisandro Dalcin       if (value[i] == '-') {
2938cc73adaaSBarry Smith         PetscCheck(i != (int)len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value);
29392d747510SLisandro Dalcin         value[i] = 0;
29402d747510SLisandro Dalcin 
29419566063dSJacob Faibussowitsch         PetscCall(PetscOptionsStringToInt(value, &start));
29422d747510SLisandro Dalcin         inc = 1;
29432d747510SLisandro Dalcin         j   = i + 1;
29442d747510SLisandro Dalcin         for (; j < (int)len; j++) {
29452d747510SLisandro Dalcin           if (value[j] == ':') {
29462d747510SLisandro Dalcin             value[j] = 0;
29472d747510SLisandro Dalcin 
29489566063dSJacob Faibussowitsch             PetscCall(PetscOptionsStringToInt(value + j + 1, &inc));
294908401ef6SPierre 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);
29502d747510SLisandro Dalcin             break;
29512d747510SLisandro Dalcin           }
29522d747510SLisandro Dalcin         }
29539566063dSJacob Faibussowitsch         PetscCall(PetscOptionsStringToInt(value + i + 1, &end));
295408401ef6SPierre 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);
29552d747510SLisandro Dalcin         nvalues = (end - start) / inc + (end - start) % inc;
2956cc73adaaSBarry 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);
29572d747510SLisandro Dalcin         for (; start < end; start += inc) {
29589371c9d4SSatish Balay           *ivalue = start;
29599371c9d4SSatish Balay           ivalue++;
29609371c9d4SSatish Balay           n++;
29612d747510SLisandro Dalcin         }
29622d747510SLisandro Dalcin         foundrange = PETSC_TRUE;
29632d747510SLisandro Dalcin         break;
29642d747510SLisandro Dalcin       }
29652d747510SLisandro Dalcin     }
29662d747510SLisandro Dalcin     if (!foundrange) {
29679566063dSJacob Faibussowitsch       PetscCall(PetscOptionsStringToInt(value, ivalue));
29682d747510SLisandro Dalcin       ivalue++;
29692d747510SLisandro Dalcin       n++;
29702d747510SLisandro Dalcin     }
29719566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
29722d747510SLisandro Dalcin   }
29739566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
29742d747510SLisandro Dalcin   *nmax = n;
29753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29762d747510SLisandro Dalcin }
29772d747510SLisandro Dalcin 
29782d747510SLisandro Dalcin /*@C
29792d747510SLisandro Dalcin   PetscOptionsGetRealArray - Gets an array of double precision values for a
2980f1a722f8SMatthew G. Knepley   particular option in the database.  The values must be separated with commas with no intervening spaces.
29812d747510SLisandro Dalcin 
29822d747510SLisandro Dalcin   Not Collective
29832d747510SLisandro Dalcin 
29842d747510SLisandro Dalcin   Input Parameters:
298520f4b53cSBarry Smith + options - options database, use `NULL` for default global database
298620f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
29876b867d5aSJose E. Roman - name    - the option one is seeking
29886b867d5aSJose E. Roman 
29892d747510SLisandro Dalcin   Output Parameters:
29902d747510SLisandro Dalcin + dvalue - the double values to return
2991f1a722f8SMatthew G. Knepley . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2992811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
29932d747510SLisandro Dalcin 
29942d747510SLisandro Dalcin   Level: beginner
29952d747510SLisandro Dalcin 
2996db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2997db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
2998db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2999c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3000db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3001db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
30022d747510SLisandro Dalcin @*/
3003d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetRealArray(PetscOptions options, const char pre[], const char name[], PetscReal dvalue[], PetscInt *nmax, PetscBool *set)
3004d71ae5a4SJacob Faibussowitsch {
30052d747510SLisandro Dalcin   const char *svalue;
30062d747510SLisandro Dalcin   char       *value;
30072d747510SLisandro Dalcin   PetscInt    n = 0;
30082d747510SLisandro Dalcin   PetscBool   flag;
30092d747510SLisandro Dalcin   PetscToken  token;
30102d747510SLisandro Dalcin 
30112d747510SLisandro Dalcin   PetscFunctionBegin;
30124f572ea9SToby Isaac   PetscAssertPointer(name, 3);
30134f572ea9SToby Isaac   PetscAssertPointer(dvalue, 4);
30144f572ea9SToby Isaac   PetscAssertPointer(nmax, 5);
30152d747510SLisandro Dalcin 
30169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
30179371c9d4SSatish Balay   if (!flag || !svalue) {
30189371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
30199371c9d4SSatish Balay     *nmax = 0;
30203ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
30219371c9d4SSatish Balay   }
30222d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
30239566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
30249566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
30252d747510SLisandro Dalcin   while (value && n < *nmax) {
30269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToReal(value, dvalue++));
30279566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
30282d747510SLisandro Dalcin     n++;
30292d747510SLisandro Dalcin   }
30309566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
30312d747510SLisandro Dalcin   *nmax = n;
30323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30332d747510SLisandro Dalcin }
30342d747510SLisandro Dalcin 
30352d747510SLisandro Dalcin /*@C
30362d747510SLisandro Dalcin   PetscOptionsGetScalarArray - Gets an array of scalars for a
3037f1a722f8SMatthew G. Knepley   particular option in the database.  The values must be separated with commas with no intervening spaces.
30382d747510SLisandro Dalcin 
30392d747510SLisandro Dalcin   Not Collective
30402d747510SLisandro Dalcin 
30412d747510SLisandro Dalcin   Input Parameters:
304220f4b53cSBarry Smith + options - options database, use `NULL` for default global database
304320f4b53cSBarry Smith . pre     - string to prepend to each name or `NULL`
30446b867d5aSJose E. Roman - name    - the option one is seeking
30456b867d5aSJose E. Roman 
30462d747510SLisandro Dalcin   Output Parameters:
30472d747510SLisandro Dalcin + dvalue - the scalar values to return
3048f1a722f8SMatthew G. Knepley . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3049811af0c4SBarry Smith - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`
30502d747510SLisandro Dalcin 
30512d747510SLisandro Dalcin   Level: beginner
30522d747510SLisandro Dalcin 
3053db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3054db781477SPatrick Sanan           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3055db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3056c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3057db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3058db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
30592d747510SLisandro Dalcin @*/
3060d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options, const char pre[], const char name[], PetscScalar dvalue[], PetscInt *nmax, PetscBool *set)
3061d71ae5a4SJacob Faibussowitsch {
30622d747510SLisandro Dalcin   const char *svalue;
30632d747510SLisandro Dalcin   char       *value;
30642d747510SLisandro Dalcin   PetscInt    n = 0;
30652d747510SLisandro Dalcin   PetscBool   flag;
30662d747510SLisandro Dalcin   PetscToken  token;
30672d747510SLisandro Dalcin 
30682d747510SLisandro Dalcin   PetscFunctionBegin;
30694f572ea9SToby Isaac   PetscAssertPointer(name, 3);
30704f572ea9SToby Isaac   PetscAssertPointer(dvalue, 4);
30714f572ea9SToby Isaac   PetscAssertPointer(nmax, 5);
30722d747510SLisandro Dalcin 
30739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
30749371c9d4SSatish Balay   if (!flag || !svalue) {
30759371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
30769371c9d4SSatish Balay     *nmax = 0;
30773ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
30789371c9d4SSatish Balay   }
30792d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
30809566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
30819566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
30822d747510SLisandro Dalcin   while (value && n < *nmax) {
30839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsStringToScalar(value, dvalue++));
30849566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
30852d747510SLisandro Dalcin     n++;
30862d747510SLisandro Dalcin   }
30879566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
30882d747510SLisandro Dalcin   *nmax = n;
30893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30902d747510SLisandro Dalcin }
309114ce751eSBarry Smith 
3092e5c89e4eSSatish Balay /*@C
3093e5c89e4eSSatish Balay   PetscOptionsGetStringArray - Gets an array of string values for a particular
3094f1a722f8SMatthew G. Knepley   option in the database. The values must be separated with commas with no intervening spaces.
3095e5c89e4eSSatish Balay 
3096cf53795eSBarry Smith   Not Collective; No Fortran Support
3097e5c89e4eSSatish Balay 
3098e5c89e4eSSatish Balay   Input Parameters:
309920f4b53cSBarry Smith + options - options database, use `NULL` for default global database
310020f4b53cSBarry Smith . pre     - string to prepend to name or `NULL`
31016b867d5aSJose E. Roman - name    - the option one is seeking
31026b867d5aSJose E. Roman 
3103e7b76fa7SPatrick Sanan   Output Parameters:
3104e5c89e4eSSatish Balay + strings - location to copy strings
3105f1a722f8SMatthew G. Knepley . nmax    - On input maximum number of strings, on output the actual number of strings found
3106811af0c4SBarry Smith - set     - `PETSC_TRUE` if found, else `PETSC_FALSE`
3107e5c89e4eSSatish Balay 
3108e5c89e4eSSatish Balay   Level: beginner
3109e5c89e4eSSatish Balay 
3110e5c89e4eSSatish Balay   Notes:
31119314d9b7SBarry Smith   The `nmax` parameter is used for both input and output.
3112e7b76fa7SPatrick Sanan 
3113e5c89e4eSSatish Balay   The user should pass in an array of pointers to char, to hold all the
3114e5c89e4eSSatish Balay   strings returned by this function.
3115e5c89e4eSSatish Balay 
3116e5c89e4eSSatish Balay   The user is responsible for deallocating the strings that are
3117cf53795eSBarry Smith   returned.
3118e5c89e4eSSatish Balay 
3119db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
3120db781477SPatrick Sanan           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3121db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3122c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3123db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3124db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
3125e5c89e4eSSatish Balay @*/
3126d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetStringArray(PetscOptions options, const char pre[], const char name[], char *strings[], PetscInt *nmax, PetscBool *set)
3127d71ae5a4SJacob Faibussowitsch {
31282d747510SLisandro Dalcin   const char *svalue;
3129e5c89e4eSSatish Balay   char       *value;
31302d747510SLisandro Dalcin   PetscInt    n = 0;
3131ace3abfcSBarry Smith   PetscBool   flag;
31329c9d3cfdSBarry Smith   PetscToken  token;
3133e5c89e4eSSatish Balay 
3134e5c89e4eSSatish Balay   PetscFunctionBegin;
31354f572ea9SToby Isaac   PetscAssertPointer(name, 3);
31364f572ea9SToby Isaac   PetscAssertPointer(strings, 4);
31374f572ea9SToby Isaac   PetscAssertPointer(nmax, 5);
3138e5c89e4eSSatish Balay 
31399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
31409371c9d4SSatish Balay   if (!flag || !svalue) {
31419371c9d4SSatish Balay     if (set) *set = PETSC_FALSE;
31429371c9d4SSatish Balay     *nmax = 0;
31433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
31449371c9d4SSatish Balay   }
31452d747510SLisandro Dalcin   if (set) *set = PETSC_TRUE;
31469566063dSJacob Faibussowitsch   PetscCall(PetscTokenCreate(svalue, ',', &token));
31479566063dSJacob Faibussowitsch   PetscCall(PetscTokenFind(token, &value));
31482d747510SLisandro Dalcin   while (value && n < *nmax) {
31499566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(value, &strings[n]));
31509566063dSJacob Faibussowitsch     PetscCall(PetscTokenFind(token, &value));
3151e5c89e4eSSatish Balay     n++;
3152e5c89e4eSSatish Balay   }
31539566063dSJacob Faibussowitsch   PetscCall(PetscTokenDestroy(&token));
3154e5c89e4eSSatish Balay   *nmax = n;
31553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3156e5c89e4eSSatish Balay }
315706824ed3SPatrick Sanan 
315806824ed3SPatrick Sanan /*@C
3159aec76313SJacob Faibussowitsch   PetscOptionsDeprecated_Private - mark an option as deprecated, optionally replacing it with `newname`
316006824ed3SPatrick Sanan 
316106824ed3SPatrick Sanan   Prints a deprecation warning, unless an option is supplied to suppress.
316206824ed3SPatrick Sanan 
31631c9f3c13SBarry Smith   Logically Collective
316406824ed3SPatrick Sanan 
316506824ed3SPatrick Sanan   Input Parameters:
3166aec76313SJacob Faibussowitsch + PetscOptionsObject - string to prepend to name or `NULL`
316706824ed3SPatrick Sanan . oldname            - the old, deprecated option
316820f4b53cSBarry Smith . newname            - the new option, or `NULL` if option is purely removed
31699f3a6782SPatrick Sanan . version            - a string describing the version of first deprecation, e.g. "3.9"
317020f4b53cSBarry Smith - info               - additional information string, or `NULL`.
317106824ed3SPatrick Sanan 
3172811af0c4SBarry Smith   Options Database Key:
317306824ed3SPatrick Sanan . -options_suppress_deprecated_warnings - do not print deprecation warnings
317406824ed3SPatrick Sanan 
317520f4b53cSBarry Smith   Level: developer
317620f4b53cSBarry Smith 
317706824ed3SPatrick Sanan   Notes:
31784ead3382SBarry Smith   If `newname` is provided then the options database will automatically check the database for `oldname`.
31794ead3382SBarry Smith 
31804ead3382SBarry Smith   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
31814ead3382SBarry Smith   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
31824ead3382SBarry Smith   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.
31834ead3382SBarry Smith 
3184811af0c4SBarry Smith   Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
318535cb6cd3SPierre Jolivet   Only the process of rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or
3186811af0c4SBarry Smith   `PetscObjectOptionsBegin()` prints the information
3187b40114eaSPatrick Sanan   If newname is provided, the old option is replaced. Otherwise, it remains
3188b40114eaSPatrick Sanan   in the options database.
31899f3a6782SPatrick Sanan   If an option is not replaced, the info argument should be used to advise the user
31909f3a6782SPatrick Sanan   on how to proceed.
31919f3a6782SPatrick Sanan   There is a limit on the length of the warning printed, so very long strings
31929f3a6782SPatrick Sanan   provided as info may be truncated.
319306824ed3SPatrick Sanan 
3194db781477SPatrick Sanan .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
319506824ed3SPatrick Sanan @*/
3196d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject, const char oldname[], const char newname[], const char version[], const char info[])
3197d71ae5a4SJacob Faibussowitsch {
319806824ed3SPatrick Sanan   PetscBool         found, quiet;
319906824ed3SPatrick Sanan   const char       *value;
320006824ed3SPatrick Sanan   const char *const quietopt = "-options_suppress_deprecated_warnings";
32019f3a6782SPatrick Sanan   char              msg[4096];
3202b0bdc838SStefano Zampini   char             *prefix  = NULL;
3203b0bdc838SStefano Zampini   PetscOptions      options = NULL;
3204b0bdc838SStefano Zampini   MPI_Comm          comm    = PETSC_COMM_SELF;
320506824ed3SPatrick Sanan 
320606824ed3SPatrick Sanan   PetscFunctionBegin;
32074f572ea9SToby Isaac   PetscAssertPointer(oldname, 2);
32084f572ea9SToby Isaac   PetscAssertPointer(version, 4);
3209b0bdc838SStefano Zampini   if (PetscOptionsObject) {
3210b0bdc838SStefano Zampini     prefix  = PetscOptionsObject->prefix;
3211b0bdc838SStefano Zampini     options = PetscOptionsObject->options;
3212b0bdc838SStefano Zampini     comm    = PetscOptionsObject->comm;
3213b0bdc838SStefano Zampini   }
32149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFindPair(options, prefix, oldname, &value, &found));
321506824ed3SPatrick Sanan   if (found) {
321606824ed3SPatrick Sanan     if (newname) {
32171baa6e33SBarry Smith       if (prefix) PetscCall(PetscOptionsPrefixPush(options, prefix));
32189566063dSJacob Faibussowitsch       PetscCall(PetscOptionsSetValue(options, newname, value));
32191baa6e33SBarry Smith       if (prefix) PetscCall(PetscOptionsPrefixPop(options));
32209566063dSJacob Faibussowitsch       PetscCall(PetscOptionsClearValue(options, oldname));
3221b40114eaSPatrick Sanan     }
322206824ed3SPatrick Sanan     quiet = PETSC_FALSE;
32239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(options, NULL, quietopt, &quiet, NULL));
322406824ed3SPatrick Sanan     if (!quiet) {
3225c6a7a370SJeremy L Thompson       PetscCall(PetscStrncpy(msg, "** PETSc DEPRECATION WARNING ** : the option ", sizeof(msg)));
3226c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, oldname, sizeof(msg)));
3227c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, " is deprecated as of version ", sizeof(msg)));
3228c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, version, sizeof(msg)));
32294bd3d7f8SBarry Smith       PetscCall(PetscStrlcat(msg, " and will be removed in a future release.\n", sizeof(msg)));
323006824ed3SPatrick Sanan       if (newname) {
32314bd3d7f8SBarry Smith         PetscCall(PetscStrlcat(msg, "   Use the option ", sizeof(msg)));
3232c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, newname, sizeof(msg)));
3233c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, " instead.", sizeof(msg)));
323406824ed3SPatrick Sanan       }
32359f3a6782SPatrick Sanan       if (info) {
3236c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, " ", sizeof(msg)));
3237c6a7a370SJeremy L Thompson         PetscCall(PetscStrlcat(msg, info, sizeof(msg)));
32389f3a6782SPatrick Sanan       }
3239c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, " (Silence this warning with ", sizeof(msg)));
3240c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, quietopt, sizeof(msg)));
3241c6a7a370SJeremy L Thompson       PetscCall(PetscStrlcat(msg, ")\n", sizeof(msg)));
32429566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%s", msg));
324306824ed3SPatrick Sanan     }
324406824ed3SPatrick Sanan   }
32453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324606824ed3SPatrick Sanan }
3247