xref: /petsc/src/sys/objects/options.c (revision fff3f5fae336d1caaec33911beb15ba4efbb662b)
1 /* Define Feature test macros to make sure atoll is available (SVr4, POSIX.1-2001, 4.3BSD, C99), not in (C89 and POSIX.1-1996) */
2 #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for atoll() */
3 
4 /*
5    These routines simplify the use of command line, file options, etc., and are used to manipulate the options database.
6    This provides the low-level interface, the high level interface is in aoptions.c
7 
8    Some routines use regular malloc and free because it cannot know  what malloc is requested with the
9    options database until it has already processed the input.
10 */
11 
12 #include <petsc/private/petscimpl.h> /*I  "petscsys.h"   I*/
13 #include <petscviewer.h>
14 #include <ctype.h>
15 #if defined(PETSC_HAVE_MALLOC_H)
16   #include <malloc.h>
17 #endif
18 #if defined(PETSC_HAVE_STRINGS_H)
19   #include <strings.h> /* strcasecmp */
20 #endif
21 
22 #if defined(PETSC_HAVE_STRCASECMP)
23   #define PetscOptNameCmp(a, b) strcasecmp(a, b)
24 #elif defined(PETSC_HAVE_STRICMP)
25   #define PetscOptNameCmp(a, b) stricmp(a, b)
26 #else
27   #define PetscOptNameCmp(a, b) Error_strcasecmp_not_found
28 #endif
29 
30 #include <petsc/private/hashtable.h>
31 
32 /* This assumes ASCII encoding and ignores locale settings */
33 /* Using tolower() is about 2X slower in microbenchmarks   */
34 static inline int PetscToLower(int c)
35 {
36   return ((c >= 'A') & (c <= 'Z')) ? c + 'a' - 'A' : c;
37 }
38 
39 /* Bob Jenkins's one at a time hash function (case-insensitive) */
40 static inline unsigned int PetscOptHash(const char key[])
41 {
42   unsigned int hash = 0;
43   while (*key) {
44     hash += PetscToLower(*key++);
45     hash += hash << 10;
46     hash ^= hash >> 6;
47   }
48   hash += hash << 3;
49   hash ^= hash >> 11;
50   hash += hash << 15;
51   return hash;
52 }
53 
54 static inline int PetscOptEqual(const char a[], const char b[])
55 {
56   return !PetscOptNameCmp(a, b);
57 }
58 
59 KHASH_INIT(HO, kh_cstr_t, int, 1, PetscOptHash, PetscOptEqual)
60 
61 #define MAXPREFIXES        25
62 #define MAXOPTIONSMONITORS 5
63 
64 const char *PetscOptionSources[] = {"code", "command line", "file", "environment"};
65 
66 // This table holds all the options set by the user
67 struct _n_PetscOptions {
68   PetscOptions previous;
69 
70   int                N;      /* number of options */
71   int                Nalloc; /* number of allocated options */
72   char             **names;  /* option names */
73   char             **values; /* option values */
74   PetscBool         *used;   /* flag option use */
75   PetscOptionSource *source; /* source for option value */
76   PetscBool          precedentProcessed;
77 
78   /* Hash table */
79   khash_t(HO) *ht;
80 
81   /* Prefixes */
82   int  prefixind;
83   int  prefixstack[MAXPREFIXES];
84   char prefix[PETSC_MAX_OPTION_NAME];
85 
86   /* Aliases */
87   int    Na;       /* number or aliases */
88   int    Naalloc;  /* number of allocated aliases */
89   char **aliases1; /* aliased */
90   char **aliases2; /* aliasee */
91 
92   /* Help */
93   PetscBool help;       /* flag whether "-help" is in the database */
94   PetscBool help_intro; /* flag whether "-help intro" is in the database */
95 
96   /* Monitors */
97   PetscBool monitorFromOptions, monitorCancel;
98   PetscErrorCode (*monitor[MAXOPTIONSMONITORS])(const char[], const char[], PetscOptionSource, void *); /* returns control to user after */
99   PetscErrorCode (*monitordestroy[MAXOPTIONSMONITORS])(void **);                                        /* callback for monitor destruction */
100   void    *monitorcontext[MAXOPTIONSMONITORS];                                                          /* to pass arbitrary user data into monitor */
101   PetscInt numbermonitors;                                                                              /* to, for instance, detect options being set */
102 };
103 
104 static PetscOptions defaultoptions = NULL; /* the options database routines query this object for options */
105 
106 /* list of options which preceed others, i.e., are processed in PetscOptionsProcessPrecedentFlags() */
107 /* these options can only take boolean values, the code will crash if given a non-boolean value */
108 static const char *precedentOptions[] = {"-petsc_ci", "-options_monitor", "-options_monitor_cancel", "-help", "-skip_petscrc"};
109 enum PetscPrecedentOption {
110   PO_CI_ENABLE,
111   PO_OPTIONS_MONITOR,
112   PO_OPTIONS_MONITOR_CANCEL,
113   PO_HELP,
114   PO_SKIP_PETSCRC,
115   PO_NUM
116 };
117 
118 PETSC_INTERN PetscErrorCode PetscOptionsSetValue_Private(PetscOptions, const char[], const char[], int *, PetscOptionSource);
119 PETSC_INTERN PetscErrorCode PetscOptionsInsertStringYAML_Private(PetscOptions, const char[], PetscOptionSource);
120 
121 /*
122     Options events monitor
123 */
124 static PetscErrorCode PetscOptionsMonitor(PetscOptions options, const char name[], const char value[], PetscOptionSource source)
125 {
126   PetscFunctionBegin;
127   if (!value) value = "";
128   if (options->monitorFromOptions) PetscCall(PetscOptionsMonitorDefault(name, value, source, NULL));
129   for (PetscInt i = 0; i < options->numbermonitors; i++) PetscCall((*options->monitor[i])(name, value, source, options->monitorcontext[i]));
130   PetscFunctionReturn(PETSC_SUCCESS);
131 }
132 
133 /*@
134    PetscOptionsCreate - Creates an empty options database.
135 
136    Logically collective
137 
138    Output Parameter:
139 .  options - Options database object
140 
141    Level: advanced
142 
143    Note:
144    Though PETSc has a concept of multiple options database the current code uses a single default `PetscOptions` object
145 
146    Developer Notes:
147    We may want eventually to pass a `MPI_Comm` to determine the ownership of the object
148 
149    This object never got developed after being introduced, it is not clear that supporting multiple `PetscOptions` objects is useful
150 
151 .seealso: `PetscOptionsDestroy()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
152 @*/
153 PetscErrorCode PetscOptionsCreate(PetscOptions *options)
154 {
155   PetscFunctionBegin;
156   PetscValidPointer(options, 1);
157   *options = (PetscOptions)calloc(1, sizeof(**options));
158   PetscCheck(*options, PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate the options database");
159   PetscFunctionReturn(PETSC_SUCCESS);
160 }
161 
162 /*@
163     PetscOptionsDestroy - Destroys an option database.
164 
165     Logically collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
166 
167   Input Parameter:
168 .  options - the `PetscOptions` object
169 
170    Level: advanced
171 
172 .seealso: `PetscOptionsInsert()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
173 @*/
174 PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
175 {
176   PetscFunctionBegin;
177   if (!*options) PetscFunctionReturn(PETSC_SUCCESS);
178   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()");
179   PetscCall(PetscOptionsClear(*options));
180   /* XXX what about monitors ? */
181   free(*options);
182   *options = NULL;
183   PetscFunctionReturn(PETSC_SUCCESS);
184 }
185 
186 /*
187     PetscOptionsCreateDefault - Creates the default global options database
188 */
189 PetscErrorCode PetscOptionsCreateDefault(void)
190 {
191   PetscFunctionBegin;
192   if (PetscUnlikely(!defaultoptions)) PetscCall(PetscOptionsCreate(&defaultoptions));
193   PetscFunctionReturn(PETSC_SUCCESS);
194 }
195 
196 /*@
197       PetscOptionsPush - Push a new `PetscOptions` object as the default provider of options
198                          Allows using different parts of a code to use different options databases
199 
200   Logically Collective
201 
202   Input Parameter:
203 .   opt - the options obtained with `PetscOptionsCreate()`
204 
205   Notes:
206   Use `PetscOptionsPop()` to return to the previous default options database
207 
208   The collectivity of this routine is complex; only the MPI ranks that call this routine will
209   have the affect of these options. If some processes that create objects call this routine and others do
210   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
211   on different ranks.
212 
213   Developer Note:
214   Though this functionality has been provided it has never been used in PETSc and might be removed.
215 
216    Level: advanced
217 
218 .seealso: `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
219 @*/
220 PetscErrorCode PetscOptionsPush(PetscOptions opt)
221 {
222   PetscFunctionBegin;
223   PetscCall(PetscOptionsCreateDefault());
224   opt->previous  = defaultoptions;
225   defaultoptions = opt;
226   PetscFunctionReturn(PETSC_SUCCESS);
227 }
228 
229 /*@
230       PetscOptionsPop - Pop the most recent `PetscOptionsPush()` to return to the previous default options
231 
232       Logically collective on whatever communicator was associated with the call to `PetscOptionsCreate()`
233 
234    Level: advanced
235 
236 .seealso: `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
237 @*/
238 PetscErrorCode PetscOptionsPop(void)
239 {
240   PetscOptions current = defaultoptions;
241 
242   PetscFunctionBegin;
243   PetscCheck(defaultoptions, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing default options");
244   PetscCheck(defaultoptions->previous, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscOptionsPop() called too many times");
245   defaultoptions    = defaultoptions->previous;
246   current->previous = NULL;
247   PetscFunctionReturn(PETSC_SUCCESS);
248 }
249 
250 /*
251     PetscOptionsDestroyDefault - Destroys the default global options database
252 */
253 PetscErrorCode PetscOptionsDestroyDefault(void)
254 {
255   PetscFunctionBegin;
256   if (!defaultoptions) PetscFunctionReturn(PETSC_SUCCESS);
257   /* Destroy any options that the user forgot to pop */
258   while (defaultoptions->previous) {
259     PetscOptions tmp = defaultoptions;
260 
261     PetscCall(PetscOptionsPop());
262     PetscCall(PetscOptionsDestroy(&tmp));
263   }
264   PetscCall(PetscOptionsDestroy(&defaultoptions));
265   PetscFunctionReturn(PETSC_SUCCESS);
266 }
267 
268 /*@C
269    PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.
270 
271    Not collective
272 
273    Input Parameter:
274 .  key - string to check if valid
275 
276    Output Parameter:
277 .  valid - `PETSC_TRUE` if a valid key
278 
279    Level: intermediate
280 @*/
281 PetscErrorCode PetscOptionsValidKey(const char key[], PetscBool *valid)
282 {
283   char *ptr;
284 
285   PetscFunctionBegin;
286   if (key) PetscValidCharPointer(key, 1);
287   PetscValidBoolPointer(valid, 2);
288   *valid = PETSC_FALSE;
289   if (!key) PetscFunctionReturn(PETSC_SUCCESS);
290   if (key[0] != '-') PetscFunctionReturn(PETSC_SUCCESS);
291   if (key[1] == '-') key++;
292   if (!isalpha((int)key[1])) PetscFunctionReturn(PETSC_SUCCESS);
293   (void)strtod(key, &ptr);
294   if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) PetscFunctionReturn(PETSC_SUCCESS);
295   *valid = PETSC_TRUE;
296   PetscFunctionReturn(PETSC_SUCCESS);
297 }
298 
299 PetscErrorCode PetscOptionsInsertString_Private(PetscOptions options, const char in_str[], PetscOptionSource source)
300 {
301   MPI_Comm   comm = PETSC_COMM_SELF;
302   char      *first, *second;
303   PetscToken token;
304 
305   PetscFunctionBegin;
306   PetscCall(PetscTokenCreate(in_str, ' ', &token));
307   PetscCall(PetscTokenFind(token, &first));
308   while (first) {
309     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
310     PetscCall(PetscStrcasecmp(first, "-options_file", &isfile));
311     PetscCall(PetscStrcasecmp(first, "-options_file_yaml", &isfileyaml));
312     PetscCall(PetscStrcasecmp(first, "-options_string_yaml", &isstringyaml));
313     PetscCall(PetscStrcasecmp(first, "-prefix_push", &ispush));
314     PetscCall(PetscStrcasecmp(first, "-prefix_pop", &ispop));
315     PetscCall(PetscOptionsValidKey(first, &key));
316     if (!key) {
317       PetscCall(PetscTokenFind(token, &first));
318     } else if (isfile) {
319       PetscCall(PetscTokenFind(token, &second));
320       PetscCall(PetscOptionsInsertFile(comm, options, second, PETSC_TRUE));
321       PetscCall(PetscTokenFind(token, &first));
322     } else if (isfileyaml) {
323       PetscCall(PetscTokenFind(token, &second));
324       PetscCall(PetscOptionsInsertFileYAML(comm, options, second, PETSC_TRUE));
325       PetscCall(PetscTokenFind(token, &first));
326     } else if (isstringyaml) {
327       PetscCall(PetscTokenFind(token, &second));
328       PetscCall(PetscOptionsInsertStringYAML_Private(options, second, source));
329       PetscCall(PetscTokenFind(token, &first));
330     } else if (ispush) {
331       PetscCall(PetscTokenFind(token, &second));
332       PetscCall(PetscOptionsPrefixPush(options, second));
333       PetscCall(PetscTokenFind(token, &first));
334     } else if (ispop) {
335       PetscCall(PetscOptionsPrefixPop(options));
336       PetscCall(PetscTokenFind(token, &first));
337     } else {
338       PetscCall(PetscTokenFind(token, &second));
339       PetscCall(PetscOptionsValidKey(second, &key));
340       if (!key) {
341         PetscCall(PetscOptionsSetValue_Private(options, first, second, NULL, source));
342         PetscCall(PetscTokenFind(token, &first));
343       } else {
344         PetscCall(PetscOptionsSetValue_Private(options, first, NULL, NULL, source));
345         first = second;
346       }
347     }
348   }
349   PetscCall(PetscTokenDestroy(&token));
350   PetscFunctionReturn(PETSC_SUCCESS);
351 }
352 
353 /*@C
354    PetscOptionsInsertString - Inserts options into the database from a string
355 
356    Logically Collective
357 
358    Input Parameters:
359 +  options - options object
360 -  in_str - string that contains options separated by blanks
361 
362    Level: intermediate
363 
364   The collectivity of this routine is complex; only the MPI ranks that call this routine will
365   have the affect of these options. If some processes that create objects call this routine and others do
366   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
367   on different ranks.
368 
369    Contributed by Boyana Norris
370 
371 .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
372           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
373           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
374           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
375           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
376           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsInsertFile()`
377 @*/
378 PetscErrorCode PetscOptionsInsertString(PetscOptions options, const char in_str[])
379 {
380   PetscFunctionBegin;
381   PetscCall(PetscOptionsInsertString_Private(options, in_str, PETSC_OPT_CODE));
382   PetscFunctionReturn(PETSC_SUCCESS);
383 }
384 
385 /*
386     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
387 */
388 static char *Petscgetline(FILE *f)
389 {
390   size_t size = 0;
391   size_t len  = 0;
392   size_t last = 0;
393   char  *buf  = NULL;
394 
395   if (feof(f)) return NULL;
396   do {
397     size += 1024;                             /* BUFSIZ is defined as "the optimal read size for this platform" */
398     buf = (char *)realloc((void *)buf, size); /* realloc(NULL,n) is the same as malloc(n) */
399     /* Actually do the read. Note that fgets puts a terminal '\0' on the
400     end of the string, so we make sure we overwrite this */
401     if (!fgets(buf + len, 1024, f)) buf[len] = 0;
402     PetscCallAbort(PETSC_COMM_SELF, PetscStrlen(buf, &len));
403     last = len - 1;
404   } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
405   if (len) return buf;
406   free(buf);
407   return NULL;
408 }
409 
410 static PetscErrorCode PetscOptionsFilename(MPI_Comm comm, const char file[], char filename[PETSC_MAX_PATH_LEN], PetscBool *yaml)
411 {
412   char fname[PETSC_MAX_PATH_LEN + 8], path[PETSC_MAX_PATH_LEN + 8], *tail;
413 
414   PetscFunctionBegin;
415   *yaml = PETSC_FALSE;
416   PetscCall(PetscStrreplace(comm, file, fname, sizeof(fname)));
417   PetscCall(PetscFixFilename(fname, path));
418   PetscCall(PetscStrendswith(path, ":yaml", yaml));
419   if (*yaml) {
420     PetscCall(PetscStrrchr(path, ':', &tail));
421     tail[-1] = 0; /* remove ":yaml" suffix from path */
422   }
423   PetscCall(PetscStrncpy(filename, path, PETSC_MAX_PATH_LEN));
424   /* check for standard YAML and JSON filename extensions */
425   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yaml", yaml));
426   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yml", yaml));
427   if (!*yaml) PetscCall(PetscStrendswith(filename, ".json", yaml));
428   if (!*yaml) { /* check file contents */
429     PetscMPIInt rank;
430     PetscCallMPI(MPI_Comm_rank(comm, &rank));
431     if (rank == 0) {
432       FILE *fh = fopen(filename, "r");
433       if (fh) {
434         char buf[6] = "";
435         if (fread(buf, 1, 6, fh) > 0) {
436           PetscCall(PetscStrncmp(buf, "%YAML ", 6, yaml));          /* check for '%YAML' tag */
437           if (!*yaml) PetscCall(PetscStrncmp(buf, "---", 3, yaml)); /* check for document start */
438         }
439         (void)fclose(fh);
440       }
441     }
442     PetscCallMPI(MPI_Bcast(yaml, 1, MPIU_BOOL, 0, comm));
443   }
444   PetscFunctionReturn(PETSC_SUCCESS);
445 }
446 
447 static PetscErrorCode PetscOptionsInsertFilePetsc(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
448 {
449   char       *string, *vstring = NULL, *astring = NULL, *packed = NULL;
450   char       *tokens[4];
451   size_t      i, len, bytes;
452   FILE       *fd;
453   PetscToken  token = NULL;
454   int         err;
455   char       *cmatch = NULL;
456   const char  cmt    = '#';
457   PetscInt    line   = 1;
458   PetscMPIInt rank, cnt = 0, acnt = 0, counts[2];
459   PetscBool   isdir, alias = PETSC_FALSE, valid;
460 
461   PetscFunctionBegin;
462   PetscCall(PetscMemzero(tokens, sizeof(tokens)));
463   PetscCallMPI(MPI_Comm_rank(comm, &rank));
464   if (rank == 0) {
465     char fpath[PETSC_MAX_PATH_LEN];
466     char fname[PETSC_MAX_PATH_LEN];
467 
468     PetscCall(PetscStrreplace(PETSC_COMM_SELF, file, fpath, sizeof(fpath)));
469     PetscCall(PetscFixFilename(fpath, fname));
470 
471     fd = fopen(fname, "r");
472     PetscCall(PetscTestDirectory(fname, 'r', &isdir));
473     PetscCheck(!isdir || !require, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified options file %s is a directory", fname);
474     if (fd && !isdir) {
475       PetscSegBuffer vseg, aseg;
476       PetscCall(PetscSegBufferCreate(1, 4000, &vseg));
477       PetscCall(PetscSegBufferCreate(1, 2000, &aseg));
478 
479       /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
480       PetscCall(PetscInfo(NULL, "Opened options file %s\n", file));
481 
482       while ((string = Petscgetline(fd))) {
483         /* eliminate comments from each line */
484         PetscCall(PetscStrchr(string, cmt, &cmatch));
485         if (cmatch) *cmatch = 0;
486         PetscCall(PetscStrlen(string, &len));
487         /* replace tabs, ^M, \n with " " */
488         for (i = 0; i < len; i++) {
489           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') string[i] = ' ';
490         }
491         PetscCall(PetscTokenCreate(string, ' ', &token));
492         PetscCall(PetscTokenFind(token, &tokens[0]));
493         if (!tokens[0]) {
494           goto destroy;
495         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
496           PetscCall(PetscTokenFind(token, &tokens[0]));
497         }
498         for (i = 1; i < 4; i++) PetscCall(PetscTokenFind(token, &tokens[i]));
499         if (!tokens[0]) {
500           goto destroy;
501         } else if (tokens[0][0] == '-') {
502           PetscCall(PetscOptionsValidKey(tokens[0], &valid));
503           PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": invalid option %s", fname, line, tokens[0]);
504           PetscCall(PetscStrlen(tokens[0], &len));
505           PetscCall(PetscSegBufferGet(vseg, len + 1, &vstring));
506           PetscCall(PetscArraycpy(vstring, tokens[0], len));
507           vstring[len] = ' ';
508           if (tokens[1]) {
509             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
510             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]);
511             PetscCall(PetscStrlen(tokens[1], &len));
512             PetscCall(PetscSegBufferGet(vseg, len + 3, &vstring));
513             vstring[0] = '"';
514             PetscCall(PetscArraycpy(vstring + 1, tokens[1], len));
515             vstring[len + 1] = '"';
516             vstring[len + 2] = ' ';
517           }
518         } else {
519           PetscCall(PetscStrcasecmp(tokens[0], "alias", &alias));
520           if (alias) {
521             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
522             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]);
523             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]);
524             PetscCall(PetscOptionsValidKey(tokens[2], &valid));
525             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]);
526             PetscCall(PetscStrlen(tokens[1], &len));
527             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
528             PetscCall(PetscArraycpy(astring, tokens[1], len));
529             astring[len] = ' ';
530 
531             PetscCall(PetscStrlen(tokens[2], &len));
532             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
533             PetscCall(PetscArraycpy(astring, tokens[2], len));
534             astring[len] = ' ';
535           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown first token in options file %s line %" PetscInt_FMT ": %s", fname, line, tokens[0]);
536         }
537         {
538           const char *extraToken = alias ? tokens[3] : tokens[2];
539           PetscCheck(!extraToken, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": extra token %s", fname, line, extraToken);
540         }
541       destroy:
542         free(string);
543         PetscCall(PetscTokenDestroy(&token));
544         alias = PETSC_FALSE;
545         line++;
546       }
547       err = fclose(fd);
548       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file %s", fname);
549       PetscCall(PetscSegBufferGetSize(aseg, &bytes)); /* size without null termination */
550       PetscCall(PetscMPIIntCast(bytes, &acnt));
551       PetscCall(PetscSegBufferGet(aseg, 1, &astring));
552       astring[0] = 0;
553       PetscCall(PetscSegBufferGetSize(vseg, &bytes)); /* size without null termination */
554       PetscCall(PetscMPIIntCast(bytes, &cnt));
555       PetscCall(PetscSegBufferGet(vseg, 1, &vstring));
556       vstring[0] = 0;
557       PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
558       PetscCall(PetscSegBufferExtractTo(aseg, packed));
559       PetscCall(PetscSegBufferExtractTo(vseg, packed + acnt + 1));
560       PetscCall(PetscSegBufferDestroy(&aseg));
561       PetscCall(PetscSegBufferDestroy(&vseg));
562     } else PetscCheck(!require, PETSC_COMM_SELF, PETSC_ERR_USER, "Unable to open options file %s", fname);
563   }
564 
565   counts[0] = acnt;
566   counts[1] = cnt;
567   err       = MPI_Bcast(counts, 2, MPI_INT, 0, comm);
568   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/");
569   acnt = counts[0];
570   cnt  = counts[1];
571   if (rank) PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
572   if (acnt || cnt) {
573     PetscCallMPI(MPI_Bcast(packed, 2 + acnt + cnt, MPI_CHAR, 0, comm));
574     astring = packed;
575     vstring = packed + acnt + 1;
576   }
577 
578   if (acnt) {
579     PetscCall(PetscTokenCreate(astring, ' ', &token));
580     PetscCall(PetscTokenFind(token, &tokens[0]));
581     while (tokens[0]) {
582       PetscCall(PetscTokenFind(token, &tokens[1]));
583       PetscCall(PetscOptionsSetAlias(options, tokens[0], tokens[1]));
584       PetscCall(PetscTokenFind(token, &tokens[0]));
585     }
586     PetscCall(PetscTokenDestroy(&token));
587   }
588 
589   if (cnt) PetscCall(PetscOptionsInsertString_Private(options, vstring, PETSC_OPT_FILE));
590   PetscCall(PetscFree(packed));
591   PetscFunctionReturn(PETSC_SUCCESS);
592 }
593 
594 /*@C
595      PetscOptionsInsertFile - Inserts options into the database from a file.
596 
597      Collective
598 
599   Input Parameters:
600 +   comm - the processes that will share the options (usually `PETSC_COMM_WORLD`)
601 .   options - options database, use NULL for default global database
602 .   file - name of file,
603            ".yml" and ".yaml" filename extensions are inserted as YAML options,
604            append ":yaml" to filename to force YAML options.
605 -   require - if `PETSC_TRUE` will generate an error if the file does not exist
606 
607   Notes:
608    Use  # for lines that are comments and which should be ignored.
609    Usually, instead of using this command, one should list the file name in the call to `PetscInitialize()`, this insures that certain options
610    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
611    calls to XXXSetFromOptions() it should not be used for options listed under PetscInitialize().
612    The collectivity of this routine is complex; only the MPI ranks in comm will
613    have the affect of these options. If some processes that create objects call this routine and others do
614    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
615    on different ranks.
616 
617   Level: developer
618 
619 .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
620           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
621           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
622           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
623           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
624           `PetscOptionsFList()`, `PetscOptionsEList()`
625 @*/
626 PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
627 {
628   char      filename[PETSC_MAX_PATH_LEN];
629   PetscBool yaml;
630 
631   PetscFunctionBegin;
632   PetscCall(PetscOptionsFilename(comm, file, filename, &yaml));
633   if (yaml) {
634     PetscCall(PetscOptionsInsertFileYAML(comm, options, filename, require));
635   } else {
636     PetscCall(PetscOptionsInsertFilePetsc(comm, options, filename, require));
637   }
638   PetscFunctionReturn(PETSC_SUCCESS);
639 }
640 
641 /*@C
642    PetscOptionsInsertArgs - Inserts options into the database from a array of strings
643 
644    Logically Collective
645 
646    Input Parameters:
647 +  options - options object
648 .  argc - the array length
649 -  args - the string array
650 
651    Level: intermediate
652 
653 .seealso: `PetscOptions`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`
654 @*/
655 PetscErrorCode PetscOptionsInsertArgs(PetscOptions options, int argc, char *args[])
656 {
657   MPI_Comm     comm  = PETSC_COMM_WORLD;
658   int          left  = PetscMax(argc, 0);
659   char *const *eargs = args;
660 
661   PetscFunctionBegin;
662   while (left) {
663     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
664     PetscCall(PetscStrcasecmp(eargs[0], "-options_file", &isfile));
665     PetscCall(PetscStrcasecmp(eargs[0], "-options_file_yaml", &isfileyaml));
666     PetscCall(PetscStrcasecmp(eargs[0], "-options_string_yaml", &isstringyaml));
667     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_push", &ispush));
668     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_pop", &ispop));
669     PetscCall(PetscOptionsValidKey(eargs[0], &key));
670     if (!key) {
671       eargs++;
672       left--;
673     } else if (isfile) {
674       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file filename option");
675       PetscCall(PetscOptionsInsertFile(comm, options, eargs[1], PETSC_TRUE));
676       eargs += 2;
677       left -= 2;
678     } else if (isfileyaml) {
679       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file_yaml filename option");
680       PetscCall(PetscOptionsInsertFileYAML(comm, options, eargs[1], PETSC_TRUE));
681       eargs += 2;
682       left -= 2;
683     } else if (isstringyaml) {
684       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing string for -options_string_yaml string option");
685       PetscCall(PetscOptionsInsertStringYAML_Private(options, eargs[1], PETSC_OPT_CODE));
686       eargs += 2;
687       left -= 2;
688     } else if (ispush) {
689       PetscCheck(left > 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option");
690       PetscCheck(eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option (prefixes cannot start with '-')");
691       PetscCall(PetscOptionsPrefixPush(options, eargs[1]));
692       eargs += 2;
693       left -= 2;
694     } else if (ispop) {
695       PetscCall(PetscOptionsPrefixPop(options));
696       eargs++;
697       left--;
698     } else {
699       PetscBool nextiskey = PETSC_FALSE;
700       if (left >= 2) PetscCall(PetscOptionsValidKey(eargs[1], &nextiskey));
701       if (left < 2 || nextiskey) {
702         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], NULL, NULL, PETSC_OPT_COMMAND_LINE));
703         eargs++;
704         left--;
705       } else {
706         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], eargs[1], NULL, PETSC_OPT_COMMAND_LINE));
707         eargs += 2;
708         left -= 2;
709       }
710     }
711   }
712   PetscFunctionReturn(PETSC_SUCCESS);
713 }
714 
715 static inline PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt, const char *val[], PetscBool set[], PetscBool *flg)
716 {
717   PetscFunctionBegin;
718   if (set[opt]) {
719     PetscCall(PetscOptionsStringToBool(val[opt], flg));
720   } else *flg = PETSC_FALSE;
721   PetscFunctionReturn(PETSC_SUCCESS);
722 }
723 
724 /* Process options with absolute precedence, these are only processed from the command line, not the environment or files */
725 static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options, int argc, char *args[], PetscBool *skip_petscrc, PetscBool *skip_petscrc_set)
726 {
727   const char *const *opt = precedentOptions;
728   const size_t       n   = PO_NUM;
729   size_t             o;
730   int                a;
731   const char       **val;
732   char             **cval;
733   PetscBool         *set, unneeded;
734 
735   PetscFunctionBegin;
736   PetscCall(PetscCalloc2(n, &cval, n, &set));
737   val = (const char **)cval;
738 
739   /* Look for options possibly set using PetscOptionsSetValue beforehand */
740   for (o = 0; o < n; o++) PetscCall(PetscOptionsFindPair(options, NULL, opt[o], &val[o], &set[o]));
741 
742   /* Loop through all args to collect last occurring value of each option */
743   for (a = 1; a < argc; a++) {
744     PetscBool valid, eq;
745 
746     PetscCall(PetscOptionsValidKey(args[a], &valid));
747     if (!valid) continue;
748     for (o = 0; o < n; o++) {
749       PetscCall(PetscStrcasecmp(args[a], opt[o], &eq));
750       if (eq) {
751         set[o] = PETSC_TRUE;
752         if (a == argc - 1 || !args[a + 1] || !args[a + 1][0] || args[a + 1][0] == '-') val[o] = NULL;
753         else val[o] = args[a + 1];
754         break;
755       }
756     }
757   }
758 
759   /* Process flags */
760   PetscCall(PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro));
761   if (options->help_intro) options->help = PETSC_TRUE;
762   else PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_HELP, val, set, &options->help));
763   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_CI_ENABLE, val, set, &unneeded));
764   /* need to manage PO_CI_ENABLE option before the PetscOptionsMonitor is turned on, so its setting is not monitored */
765   if (set[PO_CI_ENABLE]) PetscCall(PetscOptionsSetValue_Private(options, opt[PO_CI_ENABLE], val[PO_CI_ENABLE], &a, PETSC_OPT_COMMAND_LINE));
766   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL, val, set, &options->monitorCancel));
767   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR, val, set, &options->monitorFromOptions));
768   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC, val, set, skip_petscrc));
769   *skip_petscrc_set = set[PO_SKIP_PETSCRC];
770 
771   /* Store precedent options in database and mark them as used */
772   for (o = 1; o < n; o++) {
773     if (set[o]) {
774       PetscCall(PetscOptionsSetValue_Private(options, opt[o], val[o], &a, PETSC_OPT_COMMAND_LINE));
775       options->used[a] = PETSC_TRUE;
776     }
777   }
778   PetscCall(PetscFree2(cval, set));
779   options->precedentProcessed = PETSC_TRUE;
780   PetscFunctionReturn(PETSC_SUCCESS);
781 }
782 
783 static inline PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options, const char name[], PetscBool *flg)
784 {
785   PetscFunctionBegin;
786   PetscValidBoolPointer(flg, 3);
787   *flg = PETSC_FALSE;
788   if (options->precedentProcessed) {
789     for (int i = 0; i < PO_NUM; ++i) {
790       if (!PetscOptNameCmp(precedentOptions[i], name)) {
791         /* check if precedent option has been set already */
792         PetscCall(PetscOptionsFindPair(options, NULL, name, NULL, flg));
793         if (*flg) break;
794       }
795     }
796   }
797   PetscFunctionReturn(PETSC_SUCCESS);
798 }
799 
800 /*@C
801    PetscOptionsInsert - Inserts into the options database from the command line,
802                         the environmental variable and a file.
803 
804    Collective on `PETSC_COMM_WORLD`
805 
806    Input Parameters:
807 +  options - options database or NULL for the default global database
808 .  argc - count of number of command line arguments
809 .  args - the command line arguments
810 -  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
811           Use NULL or empty string to not check for code specific file.
812           Also checks ~/.petscrc, .petscrc and petscrc.
813           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
814 
815    Options Database Keys:
816 +   -options_file <filename> - read options from a file
817 -   -options_file_yaml <filename> - read options from a YAML file
818 
819    Notes:
820    Since PetscOptionsInsert() is automatically called by `PetscInitialize()`,
821    the user does not typically need to call this routine. `PetscOptionsInsert()`
822    can be called several times, adding additional entries into the database.
823 
824    See `PetscInitialize()` for options related to option database monitoring.
825 
826    Level: advanced
827 
828 .seealso: `PetscOptionsDestroy()`, `PetscOptionsView()`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`,
829           `PetscInitialize()`
830 @*/
831 PetscErrorCode PetscOptionsInsert(PetscOptions options, int *argc, char ***args, const char file[])
832 {
833   MPI_Comm    comm = PETSC_COMM_WORLD;
834   PetscMPIInt rank;
835   PetscBool   hasArgs     = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
836   PetscBool   skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;
837 
838   PetscFunctionBegin;
839   PetscCheck(!hasArgs || (args && *args), comm, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
840   PetscCallMPI(MPI_Comm_rank(comm, &rank));
841 
842   if (!options) {
843     PetscCall(PetscOptionsCreateDefault());
844     options = defaultoptions;
845   }
846   if (hasArgs) {
847     /* process options with absolute precedence */
848     PetscCall(PetscOptionsProcessPrecedentFlags(options, *argc, *args, &skipPetscrc, &skipPetscrcSet));
849     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci", &PetscCIEnabled, NULL));
850   }
851   if (file && file[0]) {
852     PetscCall(PetscOptionsInsertFile(comm, options, file, PETSC_TRUE));
853     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
854     if (!skipPetscrcSet) PetscCall(PetscOptionsGetBool(options, NULL, "-skip_petscrc", &skipPetscrc, NULL));
855   }
856   if (!skipPetscrc) {
857     char filename[PETSC_MAX_PATH_LEN];
858     PetscCall(PetscGetHomeDirectory(filename, sizeof(filename)));
859     PetscCallMPI(MPI_Bcast(filename, (int)sizeof(filename), MPI_CHAR, 0, comm));
860     if (filename[0]) PetscCall(PetscStrlcat(filename, "/.petscrc", sizeof(filename)));
861     PetscCall(PetscOptionsInsertFile(comm, options, filename, PETSC_FALSE));
862     PetscCall(PetscOptionsInsertFile(comm, options, ".petscrc", PETSC_FALSE));
863     PetscCall(PetscOptionsInsertFile(comm, options, "petscrc", PETSC_FALSE));
864   }
865 
866   /* insert environment options */
867   {
868     char  *eoptions = NULL;
869     size_t len      = 0;
870     if (rank == 0) {
871       eoptions = (char *)getenv("PETSC_OPTIONS");
872       PetscCall(PetscStrlen(eoptions, &len));
873     }
874     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
875     if (len) {
876       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
877       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
878       if (rank) eoptions[len] = 0;
879       PetscCall(PetscOptionsInsertString_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
880       if (rank) PetscCall(PetscFree(eoptions));
881     }
882   }
883 
884   /* insert YAML environment options */
885   {
886     char  *eoptions = NULL;
887     size_t len      = 0;
888     if (rank == 0) {
889       eoptions = (char *)getenv("PETSC_OPTIONS_YAML");
890       PetscCall(PetscStrlen(eoptions, &len));
891     }
892     PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
893     if (len) {
894       if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
895       PetscCallMPI(MPI_Bcast(eoptions, len, MPI_CHAR, 0, comm));
896       if (rank) eoptions[len] = 0;
897       PetscCall(PetscOptionsInsertStringYAML_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
898       if (rank) PetscCall(PetscFree(eoptions));
899     }
900   }
901 
902   /* insert command line options here because they take precedence over arguments in petscrc/environment */
903   if (hasArgs) PetscCall(PetscOptionsInsertArgs(options, *argc - 1, *args + 1));
904   PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci_portable_error_output", &PetscCIEnabledPortableErrorOutput, NULL));
905   PetscFunctionReturn(PETSC_SUCCESS);
906 }
907 
908 /* These options are not printed with PetscOptionsView() or PetscOptionsMonitor() when PetscCIEnabled is on */
909 /* TODO: get the list from the test harness, do not have it hardwired here. Maybe from gmakegentest.py */
910 static const char *PetscCIOptions[] = {
911   "malloc_debug",
912   "malloc_dump",
913   "malloc_test",
914   "malloc",
915   "nox",
916   "nox_warning",
917   "display",
918   "saws_port_auto_select",
919   "saws_port_auto_select_silent",
920   "vecscatter_mpi1",
921   "check_pointer_intensity",
922   "cuda_initialize",
923   "error_output_stdout",
924   "use_gpu_aware_mpi",
925   "checkfunctionlist",
926   "petsc_ci",
927   "petsc_ci_portable_error_output",
928 };
929 
930 static PetscBool PetscCIOption(const char *name)
931 {
932   PetscInt  idx;
933   PetscBool found;
934 
935   if (!PetscCIEnabled) return PETSC_FALSE;
936   PetscCallAbort(PETSC_COMM_SELF, PetscEListFind(PETSC_STATIC_ARRAY_LENGTH(PetscCIOptions), PetscCIOptions, name, &idx, &found));
937   return found;
938 }
939 
940 /*@C
941    PetscOptionsView - Prints the options that have been loaded. This is
942    useful for debugging purposes.
943 
944    Logically Collective
945 
946    Input Parameters:
947 +  options - options database, use NULL for default global database
948 -  viewer - must be an `PETSCVIEWERASCII` viewer
949 
950    Options Database Key:
951 .  -options_view - Activates `PetscOptionsView()` within `PetscFinalize()`
952 
953    Note:
954    Only the rank zero process of the `MPI_Comm` used to create view prints the option values. Other processes
955    may have different values but they are not printed.
956 
957    Level: advanced
958 
959 .seealso: `PetscOptionsAllUsed()`
960 @*/
961 PetscErrorCode PetscOptionsView(PetscOptions options, PetscViewer viewer)
962 {
963   PetscInt  i, N = 0;
964   PetscBool isascii;
965 
966   PetscFunctionBegin;
967   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
968   options = options ? options : defaultoptions;
969   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
970   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
971   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Only supports ASCII viewer");
972 
973   for (i = 0; i < options->N; i++) {
974     if (PetscCIOption(options->names[i])) continue;
975     N++;
976   }
977 
978   if (!N) {
979     PetscCall(PetscViewerASCIIPrintf(viewer, "#No PETSc Option Table entries\n"));
980     PetscFunctionReturn(PETSC_SUCCESS);
981   }
982 
983   PetscCall(PetscViewerASCIIPrintf(viewer, "#PETSc Option Table entries:\n"));
984   for (i = 0; i < options->N; i++) {
985     if (PetscCIOption(options->names[i])) continue;
986     if (options->values[i]) {
987       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s %s", options->names[i], options->values[i]));
988     } else {
989       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s", options->names[i]));
990     }
991     PetscCall(PetscViewerASCIIPrintf(viewer, " # (source: %s)\n", PetscOptionSources[options->source[i]]));
992   }
993   PetscCall(PetscViewerASCIIPrintf(viewer, "#End of PETSc Option Table entries\n"));
994   PetscFunctionReturn(PETSC_SUCCESS);
995 }
996 
997 /*
998    Called by error handlers to print options used in run
999 */
1000 PetscErrorCode PetscOptionsLeftError(void)
1001 {
1002   PetscInt i, nopt = 0;
1003 
1004   for (i = 0; i < defaultoptions->N; i++) {
1005     if (!defaultoptions->used[i]) {
1006       if (PetscCIOption(defaultoptions->names[i])) continue;
1007       nopt++;
1008     }
1009   }
1010   if (nopt) {
1011     PetscCall((*PetscErrorPrintf)("WARNING! There are option(s) set that were not used! Could be the program crashed before they were used or a spelling mistake, etc!\n"));
1012     for (i = 0; i < defaultoptions->N; i++) {
1013       if (!defaultoptions->used[i]) {
1014         if (PetscCIOption(defaultoptions->names[i])) continue;
1015         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]]));
1016         else PetscCall((*PetscErrorPrintf)("  Option left: name:-%s (no value) source: %s\n", defaultoptions->names[i], PetscOptionSources[defaultoptions->source[i]]));
1017       }
1018     }
1019   }
1020   return PETSC_SUCCESS;
1021 }
1022 
1023 PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
1024 {
1025   PetscInt     i, N = 0;
1026   PetscOptions options = defaultoptions;
1027 
1028   for (i = 0; i < options->N; i++) {
1029     if (PetscCIOption(options->names[i])) continue;
1030     N++;
1031   }
1032 
1033   if (N) {
1034     PetscCall((*PetscErrorPrintf)("PETSc Option Table entries:\n"));
1035   } else {
1036     PetscCall((*PetscErrorPrintf)("No PETSc Option Table entries\n"));
1037   }
1038   for (i = 0; i < options->N; i++) {
1039     if (PetscCIOption(options->names[i])) continue;
1040     if (options->values[i]) {
1041       PetscCall((*PetscErrorPrintf)("-%s %s (source: %s)\n", options->names[i], options->values[i], PetscOptionSources[options->source[i]]));
1042     } else {
1043       PetscCall((*PetscErrorPrintf)("-%s (source: %s)\n", options->names[i], PetscOptionSources[options->source[i]]));
1044     }
1045   }
1046   return PETSC_SUCCESS;
1047 }
1048 
1049 /*@C
1050    PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.
1051 
1052    Logically Collective
1053 
1054    Input Parameters:
1055 +  options - options database, or NULL for the default global database
1056 -  prefix - The string to append to the existing prefix
1057 
1058    Options Database Keys:
1059 +   -prefix_push <some_prefix_> - push the given prefix
1060 -   -prefix_pop - pop the last prefix
1061 
1062    Notes:
1063    It is common to use this in conjunction with -options_file as in
1064 
1065 $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
1066 
1067    where the files no longer require all options to be prefixed with -system2_.
1068 
1069    The collectivity of this routine is complex; only the MPI ranks that call this routine will
1070    have the affect of these options. If some processes that create objects call this routine and others do
1071    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1072    on different ranks.
1073 
1074    Level: advanced
1075 
1076 .seealso: `PetscOptionsPrefixPop()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
1077 @*/
1078 PetscErrorCode PetscOptionsPrefixPush(PetscOptions options, const char prefix[])
1079 {
1080   size_t    n;
1081   PetscInt  start;
1082   char      key[PETSC_MAX_OPTION_NAME + 1];
1083   PetscBool valid;
1084 
1085   PetscFunctionBegin;
1086   PetscValidCharPointer(prefix, 2);
1087   options = options ? options : defaultoptions;
1088   PetscCheck(options->prefixind < MAXPREFIXES, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum depth of prefix stack %d exceeded, recompile \n src/sys/objects/options.c with larger value for MAXPREFIXES", MAXPREFIXES);
1089   key[0] = '-'; /* keys must start with '-' */
1090   PetscCall(PetscStrncpy(key + 1, prefix, sizeof(key) - 1));
1091   PetscCall(PetscOptionsValidKey(key, &valid));
1092   if (!valid && options->prefixind > 0 && isdigit((int)prefix[0])) valid = PETSC_TRUE; /* If the prefix stack is not empty, make numbers a valid prefix */
1093   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" : "");
1094   start = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
1095   PetscCall(PetscStrlen(prefix, &n));
1096   PetscCheck(n + 1 <= sizeof(options->prefix) - start, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum prefix length %zu exceeded", sizeof(options->prefix));
1097   PetscCall(PetscArraycpy(options->prefix + start, prefix, n + 1));
1098   options->prefixstack[options->prefixind++] = start + n;
1099   PetscFunctionReturn(PETSC_SUCCESS);
1100 }
1101 
1102 /*@C
1103    PetscOptionsPrefixPop - Remove the latest options prefix, see `PetscOptionsPrefixPush()` for details
1104 
1105    Logically Collective on the `MPI_Comm` used when called `PetscOptionsPrefixPush()`
1106 
1107   Input Parameter:
1108 .  options - options database, or NULL for the default global database
1109 
1110    Level: advanced
1111 
1112 .seealso: `PetscOptionsPrefixPush()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
1113 @*/
1114 PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1115 {
1116   PetscInt offset;
1117 
1118   PetscFunctionBegin;
1119   options = options ? options : defaultoptions;
1120   PetscCheck(options->prefixind >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "More prefixes popped than pushed");
1121   options->prefixind--;
1122   offset                  = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
1123   options->prefix[offset] = 0;
1124   PetscFunctionReturn(PETSC_SUCCESS);
1125 }
1126 
1127 /*@C
1128     PetscOptionsClear - Removes all options form the database leaving it empty.
1129 
1130     Logically Collective
1131 
1132   Input Parameter:
1133 .  options - options database, use NULL for the default global database
1134 
1135    The collectivity of this routine is complex; only the MPI ranks that call this routine will
1136    have the affect of these options. If some processes that create objects call this routine and others do
1137    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1138    on different ranks.
1139 
1140    Level: developer
1141 
1142 .seealso: `PetscOptionsInsert()`
1143 @*/
1144 PetscErrorCode PetscOptionsClear(PetscOptions options)
1145 {
1146   PetscInt i;
1147 
1148   PetscFunctionBegin;
1149   options = options ? options : defaultoptions;
1150   if (!options) PetscFunctionReturn(PETSC_SUCCESS);
1151 
1152   for (i = 0; i < options->N; i++) {
1153     if (options->names[i]) free(options->names[i]);
1154     if (options->values[i]) free(options->values[i]);
1155   }
1156   options->N = 0;
1157   free(options->names);
1158   free(options->values);
1159   free(options->used);
1160   free(options->source);
1161   options->names  = NULL;
1162   options->values = NULL;
1163   options->used   = NULL;
1164   options->source = NULL;
1165   options->Nalloc = 0;
1166 
1167   for (i = 0; i < options->Na; i++) {
1168     free(options->aliases1[i]);
1169     free(options->aliases2[i]);
1170   }
1171   options->Na = 0;
1172   free(options->aliases1);
1173   free(options->aliases2);
1174   options->aliases1 = options->aliases2 = NULL;
1175   options->Naalloc                      = 0;
1176 
1177   /* destroy hash table */
1178   kh_destroy(HO, options->ht);
1179   options->ht = NULL;
1180 
1181   options->prefixind  = 0;
1182   options->prefix[0]  = 0;
1183   options->help       = PETSC_FALSE;
1184   options->help_intro = PETSC_FALSE;
1185   PetscFunctionReturn(PETSC_SUCCESS);
1186 }
1187 
1188 /*@C
1189    PetscOptionsSetAlias - Makes a key and alias for another key
1190 
1191    Logically Collective
1192 
1193    Input Parameters:
1194 +  options - options database, or NULL for default global database
1195 .  newname - the alias
1196 -  oldname - the name that alias will refer to
1197 
1198    Level: advanced
1199 
1200    The collectivity of this routine is complex; only the MPI ranks that call this routine will
1201    have the affect of these options. If some processes that create objects call this routine and others do
1202    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1203    on different ranks.
1204 
1205 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1206           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1207           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1208           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1209           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1210           `PetscOptionsFList()`, `PetscOptionsEList()`
1211 @*/
1212 PetscErrorCode PetscOptionsSetAlias(PetscOptions options, const char newname[], const char oldname[])
1213 {
1214   size_t    len;
1215   PetscBool valid;
1216 
1217   PetscFunctionBegin;
1218   PetscValidCharPointer(newname, 2);
1219   PetscValidCharPointer(oldname, 3);
1220   options = options ? options : defaultoptions;
1221   PetscCall(PetscOptionsValidKey(newname, &valid));
1222   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliased option %s", newname);
1223   PetscCall(PetscOptionsValidKey(oldname, &valid));
1224   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliasee option %s", oldname);
1225 
1226   if (options->Na == options->Naalloc) {
1227     char **tmpA1, **tmpA2;
1228 
1229     options->Naalloc = PetscMax(4, options->Naalloc * 2);
1230     tmpA1            = (char **)malloc(options->Naalloc * sizeof(char *));
1231     tmpA2            = (char **)malloc(options->Naalloc * sizeof(char *));
1232     for (int i = 0; i < options->Na; ++i) {
1233       tmpA1[i] = options->aliases1[i];
1234       tmpA2[i] = options->aliases2[i];
1235     }
1236     free(options->aliases1);
1237     free(options->aliases2);
1238     options->aliases1 = tmpA1;
1239     options->aliases2 = tmpA2;
1240   }
1241   newname++;
1242   oldname++;
1243   PetscCall(PetscStrlen(newname, &len));
1244   options->aliases1[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1245   PetscCall(PetscStrncpy(options->aliases1[options->Na], newname, len + 1));
1246   PetscCall(PetscStrlen(oldname, &len));
1247   options->aliases2[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1248   PetscCall(PetscStrncpy(options->aliases2[options->Na], oldname, len + 1));
1249   ++options->Na;
1250   PetscFunctionReturn(PETSC_SUCCESS);
1251 }
1252 
1253 /*@C
1254    PetscOptionsSetValue - Sets an option name-value pair in the options
1255    database, overriding whatever is already present.
1256 
1257    Logically Collective
1258 
1259    Input Parameters:
1260 +  options - options database, use NULL for the default global database
1261 .  name - name of option, this SHOULD have the - prepended
1262 -  value - the option value (not used for all options, so can be NULL)
1263 
1264    Level: intermediate
1265 
1266    Note:
1267    This function can be called BEFORE `PetscInitialize()`
1268 
1269    The collectivity of this routine is complex; only the MPI ranks that call this routine will
1270    have the affect of these options. If some processes that create objects call this routine and others do
1271    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1272    on different ranks.
1273 
1274    Developers Note: Uses malloc() directly because PETSc may not be initialized yet.
1275 
1276 .seealso: `PetscOptionsInsert()`, `PetscOptionsClearValue()`
1277 @*/
1278 PetscErrorCode PetscOptionsSetValue(PetscOptions options, const char name[], const char value[])
1279 {
1280   PetscFunctionBegin;
1281   PetscCall(PetscOptionsSetValue_Private(options, name, value, NULL, PETSC_OPT_CODE));
1282   PetscFunctionReturn(PETSC_SUCCESS);
1283 }
1284 
1285 PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options, const char name[], const char value[], int *pos, PetscOptionSource source)
1286 {
1287   size_t    len;
1288   int       n, i;
1289   char    **names;
1290   char      fullname[PETSC_MAX_OPTION_NAME] = "";
1291   PetscBool flg;
1292 
1293   PetscFunctionBegin;
1294   if (!options) {
1295     PetscCall(PetscOptionsCreateDefault());
1296     options = defaultoptions;
1297   }
1298   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "name %s must start with '-'", name);
1299 
1300   PetscCall(PetscOptionsSkipPrecedent(options, name, &flg));
1301   if (flg) PetscFunctionReturn(PETSC_SUCCESS);
1302 
1303   name++; /* skip starting dash */
1304 
1305   if (options->prefixind > 0) {
1306     strncpy(fullname, options->prefix, sizeof(fullname));
1307     fullname[sizeof(fullname) - 1] = 0;
1308     strncat(fullname, name, sizeof(fullname) - strlen(fullname) - 1);
1309     fullname[sizeof(fullname) - 1] = 0;
1310     name                           = fullname;
1311   }
1312 
1313   /* check against aliases */
1314   for (i = 0; i < options->Na; i++) {
1315     int result = PetscOptNameCmp(options->aliases1[i], name);
1316     if (!result) {
1317       name = options->aliases2[i];
1318       break;
1319     }
1320   }
1321 
1322   /* slow search */
1323   n     = options->N;
1324   names = options->names;
1325   for (i = 0; i < options->N; i++) {
1326     int result = PetscOptNameCmp(names[i], name);
1327     if (!result) {
1328       n = i;
1329       goto setvalue;
1330     } else if (result > 0) {
1331       n = i;
1332       break;
1333     }
1334   }
1335   if (options->N == options->Nalloc) {
1336     char             **names, **values;
1337     PetscBool         *used;
1338     PetscOptionSource *source;
1339 
1340     options->Nalloc = PetscMax(10, options->Nalloc * 2);
1341     names           = (char **)malloc(options->Nalloc * sizeof(char *));
1342     values          = (char **)malloc(options->Nalloc * sizeof(char *));
1343     used            = (PetscBool *)malloc(options->Nalloc * sizeof(PetscBool));
1344     source          = (PetscOptionSource *)malloc(options->Nalloc * sizeof(PetscOptionSource));
1345     for (int i = 0; i < options->N; ++i) {
1346       names[i]  = options->names[i];
1347       values[i] = options->values[i];
1348       used[i]   = options->used[i];
1349       source[i] = options->source[i];
1350     }
1351     free(options->names);
1352     free(options->values);
1353     free(options->used);
1354     free(options->source);
1355     options->names  = names;
1356     options->values = values;
1357     options->used   = used;
1358     options->source = source;
1359   }
1360 
1361   /* shift remaining values up 1 */
1362   for (i = options->N; i > n; i--) {
1363     options->names[i]  = options->names[i - 1];
1364     options->values[i] = options->values[i - 1];
1365     options->used[i]   = options->used[i - 1];
1366     options->source[i] = options->source[i - 1];
1367   }
1368   options->names[n]  = NULL;
1369   options->values[n] = NULL;
1370   options->used[n]   = PETSC_FALSE;
1371   options->source[n] = PETSC_OPT_CODE;
1372   options->N++;
1373 
1374   /* destroy hash table */
1375   kh_destroy(HO, options->ht);
1376   options->ht = NULL;
1377 
1378   /* set new name */
1379   len               = strlen(name);
1380   options->names[n] = (char *)malloc((len + 1) * sizeof(char));
1381   PetscCheck(options->names[n], PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate option name");
1382   strcpy(options->names[n], name);
1383 
1384 setvalue:
1385   /* set new value */
1386   if (options->values[n]) free(options->values[n]);
1387   len = value ? strlen(value) : 0;
1388   if (len) {
1389     options->values[n] = (char *)malloc((len + 1) * sizeof(char));
1390     if (!options->values[n]) return PETSC_ERR_MEM;
1391     strcpy(options->values[n], value);
1392   } else {
1393     options->values[n] = NULL;
1394   }
1395   options->source[n] = source;
1396 
1397   /* handle -help so that it can be set from anywhere */
1398   if (!PetscOptNameCmp(name, "help")) {
1399     options->help       = PETSC_TRUE;
1400     options->help_intro = (value && !PetscOptNameCmp(value, "intro")) ? PETSC_TRUE : PETSC_FALSE;
1401     options->used[n]    = PETSC_TRUE;
1402   }
1403 
1404   PetscCall(PetscOptionsMonitor(options, name, value, source));
1405   if (pos) *pos = n;
1406   PetscFunctionReturn(PETSC_SUCCESS);
1407 }
1408 
1409 /*@C
1410    PetscOptionsClearValue - Clears an option name-value pair in the options
1411    database, overriding whatever is already present.
1412 
1413    Logically Collective
1414 
1415    Input Parameters:
1416 +  options - options database, use NULL for the default global database
1417 -  name - name of option, this SHOULD have the - prepended
1418 
1419    Level: intermediate
1420 
1421    Note:
1422    The collectivity of this routine is complex; only the MPI ranks that call this routine will
1423    have the affect of these options. If some processes that create objects call this routine and others do
1424    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1425    on different ranks.
1426 
1427 .seealso: `PetscOptionsInsert()`
1428 @*/
1429 PetscErrorCode PetscOptionsClearValue(PetscOptions options, const char name[])
1430 {
1431   int    N, n, i;
1432   char **names;
1433 
1434   PetscFunctionBegin;
1435   options = options ? options : defaultoptions;
1436   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1437   if (!PetscOptNameCmp(name, "-help")) options->help = options->help_intro = PETSC_FALSE;
1438 
1439   name++; /* skip starting dash */
1440 
1441   /* slow search */
1442   N = n = options->N;
1443   names = options->names;
1444   for (i = 0; i < N; i++) {
1445     int result = PetscOptNameCmp(names[i], name);
1446     if (!result) {
1447       n = i;
1448       break;
1449     } else if (result > 0) {
1450       n = N;
1451       break;
1452     }
1453   }
1454   if (n == N) PetscFunctionReturn(PETSC_SUCCESS); /* it was not present */
1455 
1456   /* remove name and value */
1457   if (options->names[n]) free(options->names[n]);
1458   if (options->values[n]) free(options->values[n]);
1459   /* shift remaining values down 1 */
1460   for (i = n; i < N - 1; i++) {
1461     options->names[i]  = options->names[i + 1];
1462     options->values[i] = options->values[i + 1];
1463     options->used[i]   = options->used[i + 1];
1464     options->source[i] = options->source[i + 1];
1465   }
1466   options->N--;
1467 
1468   /* destroy hash table */
1469   kh_destroy(HO, options->ht);
1470   options->ht = NULL;
1471 
1472   PetscCall(PetscOptionsMonitor(options, name, NULL, PETSC_OPT_CODE));
1473   PetscFunctionReturn(PETSC_SUCCESS);
1474 }
1475 
1476 /*@C
1477    PetscOptionsFindPair - Gets an option name-value pair from the options database.
1478 
1479    Not Collective
1480 
1481    Input Parameters:
1482 +  options - options database, use NULL for the default global database
1483 .  pre - the string to prepend to the name or NULL, this SHOULD NOT have the "-" prepended
1484 -  name - name of option, this SHOULD have the "-" prepended
1485 
1486    Output Parameters:
1487 +  value - the option value (optional, not used for all options)
1488 -  set - whether the option is set (optional)
1489 
1490    Note:
1491    Each process may find different values or no value depending on how options were inserted into the database
1492 
1493    Level: developer
1494 
1495 .seealso: `PetscOptionsSetValue()`, `PetscOptionsClearValue()`
1496 @*/
1497 PetscErrorCode PetscOptionsFindPair(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1498 {
1499   char      buf[PETSC_MAX_OPTION_NAME];
1500   PetscBool usehashtable = PETSC_TRUE;
1501   PetscBool matchnumbers = PETSC_TRUE;
1502 
1503   PetscFunctionBegin;
1504   options = options ? options : defaultoptions;
1505   PetscCheck(!pre || !PetscUnlikely(pre[0] == '-'), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1506   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1507 
1508   name++; /* skip starting dash */
1509 
1510   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1511   if (pre && pre[0]) {
1512     char *ptr = buf;
1513     if (name[0] == '-') {
1514       *ptr++ = '-';
1515       name++;
1516     }
1517     PetscCall(PetscStrncpy(ptr, pre, buf + sizeof(buf) - ptr));
1518     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
1519     name = buf;
1520   }
1521 
1522   if (PetscDefined(USE_DEBUG)) {
1523     PetscBool valid;
1524     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
1525     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
1526     PetscCall(PetscOptionsValidKey(key, &valid));
1527     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1528   }
1529 
1530   if (!options->ht && usehashtable) {
1531     int          i, ret;
1532     khiter_t     it;
1533     khash_t(HO) *ht;
1534     ht = kh_init(HO);
1535     PetscCheck(ht, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1536     ret = kh_resize(HO, ht, options->N * 2); /* twice the required size to reduce risk of collisions */
1537     PetscCheck(!ret, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1538     for (i = 0; i < options->N; i++) {
1539       it = kh_put(HO, ht, options->names[i], &ret);
1540       PetscCheck(ret == 1, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1541       kh_val(ht, it) = i;
1542     }
1543     options->ht = ht;
1544   }
1545 
1546   if (usehashtable) { /* fast search */
1547     khash_t(HO) *ht = options->ht;
1548     khiter_t     it = kh_get(HO, ht, name);
1549     if (it != kh_end(ht)) {
1550       int i            = kh_val(ht, it);
1551       options->used[i] = PETSC_TRUE;
1552       if (value) *value = options->values[i];
1553       if (set) *set = PETSC_TRUE;
1554       PetscFunctionReturn(PETSC_SUCCESS);
1555     }
1556   } else { /* slow search */
1557     int i, N = options->N;
1558     for (i = 0; i < N; i++) {
1559       int result = PetscOptNameCmp(options->names[i], name);
1560       if (!result) {
1561         options->used[i] = PETSC_TRUE;
1562         if (value) *value = options->values[i];
1563         if (set) *set = PETSC_TRUE;
1564         PetscFunctionReturn(PETSC_SUCCESS);
1565       } else if (result > 0) {
1566         break;
1567       }
1568     }
1569   }
1570 
1571   /*
1572    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
1573    Maybe this special lookup mode should be enabled on request with a push/pop API.
1574    The feature of matching _%d_ used sparingly in the codebase.
1575    */
1576   if (matchnumbers) {
1577     int i, j, cnt = 0, locs[16], loce[16];
1578     /* determine the location and number of all _%d_ in the key */
1579     for (i = 0; name[i]; i++) {
1580       if (name[i] == '_') {
1581         for (j = i + 1; name[j]; j++) {
1582           if (name[j] >= '0' && name[j] <= '9') continue;
1583           if (name[j] == '_' && j > i + 1) { /* found a number */
1584             locs[cnt]   = i + 1;
1585             loce[cnt++] = j + 1;
1586           }
1587           i = j - 1;
1588           break;
1589         }
1590       }
1591     }
1592     for (i = 0; i < cnt; i++) {
1593       PetscBool found;
1594       char      opt[PETSC_MAX_OPTION_NAME + 1] = "-", tmp[PETSC_MAX_OPTION_NAME];
1595       PetscCall(PetscStrncpy(tmp, name, PetscMin((size_t)(locs[i] + 1), sizeof(tmp))));
1596       PetscCall(PetscStrlcat(opt, tmp, sizeof(opt)));
1597       PetscCall(PetscStrlcat(opt, name + loce[i], sizeof(opt)));
1598       PetscCall(PetscOptionsFindPair(options, NULL, opt, value, &found));
1599       if (found) {
1600         if (set) *set = PETSC_TRUE;
1601         PetscFunctionReturn(PETSC_SUCCESS);
1602       }
1603     }
1604   }
1605 
1606   if (set) *set = PETSC_FALSE;
1607   PetscFunctionReturn(PETSC_SUCCESS);
1608 }
1609 
1610 /* Check whether any option begins with pre+name */
1611 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1612 {
1613   char buf[PETSC_MAX_OPTION_NAME];
1614   int  numCnt = 0, locs[16], loce[16];
1615 
1616   PetscFunctionBegin;
1617   options = options ? options : defaultoptions;
1618   PetscCheck(!pre || pre[0] != '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1619   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1620 
1621   name++; /* skip starting dash */
1622 
1623   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1624   if (pre && pre[0]) {
1625     char *ptr = buf;
1626     if (name[0] == '-') {
1627       *ptr++ = '-';
1628       name++;
1629     }
1630     PetscCall(PetscStrncpy(ptr, pre, sizeof(buf) - ((ptr == buf) ? 0 : 1)));
1631     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
1632     name = buf;
1633   }
1634 
1635   if (PetscDefined(USE_DEBUG)) {
1636     PetscBool valid;
1637     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
1638     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
1639     PetscCall(PetscOptionsValidKey(key, &valid));
1640     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1641   }
1642 
1643   /* determine the location and number of all _%d_ in the key */
1644   {
1645     int i, j;
1646     for (i = 0; name[i]; i++) {
1647       if (name[i] == '_') {
1648         for (j = i + 1; name[j]; j++) {
1649           if (name[j] >= '0' && name[j] <= '9') continue;
1650           if (name[j] == '_' && j > i + 1) { /* found a number */
1651             locs[numCnt]   = i + 1;
1652             loce[numCnt++] = j + 1;
1653           }
1654           i = j - 1;
1655           break;
1656         }
1657       }
1658     }
1659   }
1660 
1661   /* slow search */
1662   for (int c = -1; c < numCnt; ++c) {
1663     char   opt[PETSC_MAX_OPTION_NAME + 2] = "";
1664     size_t len;
1665 
1666     if (c < 0) {
1667       PetscCall(PetscStrncpy(opt, name, sizeof(opt)));
1668     } else {
1669       PetscCall(PetscStrncpy(opt, name, PetscMin((size_t)(locs[c] + 1), sizeof(opt))));
1670       PetscCall(PetscStrlcat(opt, name + loce[c], sizeof(opt) - 1));
1671     }
1672     PetscCall(PetscStrlen(opt, &len));
1673     for (int i = 0; i < options->N; i++) {
1674       PetscBool match;
1675 
1676       PetscCall(PetscStrncmp(options->names[i], opt, len, &match));
1677       if (match) {
1678         options->used[i] = PETSC_TRUE;
1679         if (value) *value = options->values[i];
1680         if (set) *set = PETSC_TRUE;
1681         PetscFunctionReturn(PETSC_SUCCESS);
1682       }
1683     }
1684   }
1685 
1686   if (set) *set = PETSC_FALSE;
1687   PetscFunctionReturn(PETSC_SUCCESS);
1688 }
1689 
1690 /*@C
1691    PetscOptionsReject - Generates an error if a certain option is given.
1692 
1693    Not Collective
1694 
1695    Input Parameters:
1696 +  options - options database, use NULL for default global database
1697 .  pre - the option prefix (may be NULL)
1698 .  name - the option name one is seeking
1699 -  mess - error message (may be NULL)
1700 
1701    Level: advanced
1702 
1703 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `OptionsHasName()`,
1704           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1705           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1706           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1707           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1708           `PetscOptionsFList()`, `PetscOptionsEList()`
1709 @*/
1710 PetscErrorCode PetscOptionsReject(PetscOptions options, const char pre[], const char name[], const char mess[])
1711 {
1712   PetscBool flag = PETSC_FALSE;
1713 
1714   PetscFunctionBegin;
1715   PetscCall(PetscOptionsHasName(options, pre, name, &flag));
1716   if (flag) {
1717     PetscCheck(!mess || !mess[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s with %s", pre ? pre : "", name + 1, mess);
1718     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s", pre ? pre : "", name + 1);
1719   }
1720   PetscFunctionReturn(PETSC_SUCCESS);
1721 }
1722 
1723 /*@C
1724    PetscOptionsHasHelp - Determines whether the "-help" option is in the database.
1725 
1726    Not Collective
1727 
1728    Input Parameter:
1729 .  options - options database, use NULL for default global database
1730 
1731    Output Parameter:
1732 .  set - `PETSC_TRUE` if found else `PETSC_FALSE`.
1733 
1734    Level: advanced
1735 
1736 .seealso: `PetscOptionsHasName()`
1737 @*/
1738 PetscErrorCode PetscOptionsHasHelp(PetscOptions options, PetscBool *set)
1739 {
1740   PetscFunctionBegin;
1741   PetscValidBoolPointer(set, 2);
1742   options = options ? options : defaultoptions;
1743   *set    = options->help;
1744   PetscFunctionReturn(PETSC_SUCCESS);
1745 }
1746 
1747 PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options, PetscBool *set)
1748 {
1749   PetscFunctionBegin;
1750   PetscValidBoolPointer(set, 2);
1751   options = options ? options : defaultoptions;
1752   *set    = options->help_intro;
1753   PetscFunctionReturn(PETSC_SUCCESS);
1754 }
1755 
1756 /*@C
1757    PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or Boolean, even
1758                       if its value is set to false.
1759 
1760    Not Collective
1761 
1762    Input Parameters:
1763 +  options - options database, use NULL for default global database
1764 .  pre - string to prepend to the name or NULL
1765 -  name - the option one is seeking
1766 
1767    Output Parameter:
1768 .  set - `PETSC_TRUE` if found else `PETSC_FALSE`.
1769 
1770    Level: beginner
1771 
1772    Note:
1773    In many cases you probably want to use `PetscOptionsGetBool()` instead of calling this, to allowing toggling values.
1774 
1775 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
1776           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1777           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1778           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1779           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1780           `PetscOptionsFList()`, `PetscOptionsEList()`
1781 @*/
1782 PetscErrorCode PetscOptionsHasName(PetscOptions options, const char pre[], const char name[], PetscBool *set)
1783 {
1784   const char *value;
1785   PetscBool   flag;
1786 
1787   PetscFunctionBegin;
1788   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
1789   if (set) *set = flag;
1790   PetscFunctionReturn(PETSC_SUCCESS);
1791 }
1792 
1793 /*@C
1794    PetscOptionsGetAll - Lists all the options the program was run with in a single string.
1795 
1796    Not Collective
1797 
1798    Input Parameter:
1799 .  options - the options database, use NULL for the default global database
1800 
1801    Output Parameter:
1802 .  copts - pointer where string pointer is stored
1803 
1804    Notes:
1805     The array and each entry in the array should be freed with `PetscFree()`
1806 
1807     Each process may have different values depending on how the options were inserted into the database
1808 
1809    Level: advanced
1810 
1811 .seealso: `PetscOptionsAllUsed()`, `PetscOptionsView()`, `PetscOptionsPush()`, `PetscOptionsPop()`
1812 @*/
1813 PetscErrorCode PetscOptionsGetAll(PetscOptions options, char *copts[])
1814 {
1815   PetscInt i;
1816   size_t   len = 1, lent = 0;
1817   char    *coptions = NULL;
1818 
1819   PetscFunctionBegin;
1820   PetscValidPointer(copts, 2);
1821   options = options ? options : defaultoptions;
1822   /* count the length of the required string */
1823   for (i = 0; i < options->N; i++) {
1824     PetscCall(PetscStrlen(options->names[i], &lent));
1825     len += 2 + lent;
1826     if (options->values[i]) {
1827       PetscCall(PetscStrlen(options->values[i], &lent));
1828       len += 1 + lent;
1829     }
1830   }
1831   PetscCall(PetscMalloc1(len, &coptions));
1832   coptions[0] = 0;
1833   for (i = 0; i < options->N; i++) {
1834     PetscCall(PetscStrlcat(coptions, "-", len));
1835     PetscCall(PetscStrlcat(coptions, options->names[i], len));
1836     PetscCall(PetscStrlcat(coptions, " ", len));
1837     if (options->values[i]) {
1838       PetscCall(PetscStrlcat(coptions, options->values[i], len));
1839       PetscCall(PetscStrlcat(coptions, " ", len));
1840     }
1841   }
1842   *copts = coptions;
1843   PetscFunctionReturn(PETSC_SUCCESS);
1844 }
1845 
1846 /*@C
1847    PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database
1848 
1849    Not Collective
1850 
1851    Input Parameters:
1852 +  options - options database, use NULL for default global database
1853 -  name - string name of option
1854 
1855    Output Parameter:
1856 .  used - `PETSC_TRUE` if the option was used, otherwise false, including if option was not found in options database
1857 
1858    Level: advanced
1859 
1860    Note:
1861    The value returned may be different on each process and depends on which options have been processed
1862    on the given process
1863 
1864 .seealso: `PetscOptionsView()`, `PetscOptionsLeft()`, `PetscOptionsAllUsed()`
1865 @*/
1866 PetscErrorCode PetscOptionsUsed(PetscOptions options, const char *name, PetscBool *used)
1867 {
1868   PetscInt i;
1869 
1870   PetscFunctionBegin;
1871   PetscValidCharPointer(name, 2);
1872   PetscValidBoolPointer(used, 3);
1873   options = options ? options : defaultoptions;
1874   *used   = PETSC_FALSE;
1875   for (i = 0; i < options->N; i++) {
1876     PetscCall(PetscStrcasecmp(options->names[i], name, used));
1877     if (*used) {
1878       *used = options->used[i];
1879       break;
1880     }
1881   }
1882   PetscFunctionReturn(PETSC_SUCCESS);
1883 }
1884 
1885 /*@
1886    PetscOptionsAllUsed - Returns a count of the number of options in the
1887    database that have never been selected.
1888 
1889    Not Collective
1890 
1891    Input Parameter:
1892 .  options - options database, use NULL for default global database
1893 
1894    Output Parameter:
1895 .  N - count of options not used
1896 
1897    Level: advanced
1898 
1899    Note:
1900    The value returned may be different on each process and depends on which options have been processed
1901    on the given process
1902 
1903 .seealso: `PetscOptionsView()`
1904 @*/
1905 PetscErrorCode PetscOptionsAllUsed(PetscOptions options, PetscInt *N)
1906 {
1907   PetscInt i, n = 0;
1908 
1909   PetscFunctionBegin;
1910   PetscValidIntPointer(N, 2);
1911   options = options ? options : defaultoptions;
1912   for (i = 0; i < options->N; i++) {
1913     if (!options->used[i]) n++;
1914   }
1915   *N = n;
1916   PetscFunctionReturn(PETSC_SUCCESS);
1917 }
1918 
1919 /*@
1920    PetscOptionsLeft - Prints to screen any options that were set and never used.
1921 
1922    Not Collective
1923 
1924    Input Parameter:
1925 .  options - options database; use NULL for default global database
1926 
1927    Options Database Key:
1928 .  -options_left - activates `PetscOptionsAllUsed()` within `PetscFinalize()`
1929 
1930    Notes:
1931       This is rarely used directly, it is called by `PetscFinalize()` in debug more or if -options_left
1932       is passed otherwise to help users determine possible mistakes in their usage of options. This
1933       only prints values on process zero of `PETSC_COMM_WORLD`.
1934 
1935       Other processes depending the objects
1936       used may have different options that are left unused.
1937 
1938    Level: advanced
1939 
1940 .seealso: `PetscOptionsAllUsed()`
1941 @*/
1942 PetscErrorCode PetscOptionsLeft(PetscOptions options)
1943 {
1944   PetscInt     i;
1945   PetscInt     cnt = 0;
1946   PetscOptions toptions;
1947 
1948   PetscFunctionBegin;
1949   toptions = options ? options : defaultoptions;
1950   for (i = 0; i < toptions->N; i++) {
1951     if (!toptions->used[i]) {
1952       if (PetscCIOption(toptions->names[i])) continue;
1953       if (toptions->values[i]) {
1954         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s value: %s source: %s\n", toptions->names[i], toptions->values[i], PetscOptionSources[toptions->source[i]]));
1955       } else {
1956         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s (no value) source: %s\n", toptions->names[i], PetscOptionSources[toptions->source[i]]));
1957       }
1958     }
1959   }
1960   if (!options) {
1961     toptions = defaultoptions;
1962     while (toptions->previous) {
1963       cnt++;
1964       toptions = toptions->previous;
1965     }
1966     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));
1967   }
1968   PetscFunctionReturn(PETSC_SUCCESS);
1969 }
1970 
1971 /*@C
1972    PetscOptionsLeftGet - Returns all options that were set and never used.
1973 
1974    Not Collective
1975 
1976    Input Parameter:
1977 .  options - options database, use NULL for default global database
1978 
1979    Output Parameters:
1980 +  N - count of options not used
1981 .  names - names of options not used
1982 -  values - values of options not used
1983 
1984    Level: advanced
1985 
1986    Notes:
1987    Users should call `PetscOptionsLeftRestore()` to free the memory allocated in this routine
1988 
1989    The value returned may be different on each process and depends on which options have been processed
1990    on the given process
1991 
1992 .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`
1993 @*/
1994 PetscErrorCode PetscOptionsLeftGet(PetscOptions options, PetscInt *N, char **names[], char **values[])
1995 {
1996   PetscInt i, n;
1997 
1998   PetscFunctionBegin;
1999   if (N) PetscValidIntPointer(N, 2);
2000   if (names) PetscValidPointer(names, 3);
2001   if (values) PetscValidPointer(values, 4);
2002   options = options ? options : defaultoptions;
2003 
2004   /* The number of unused PETSc options */
2005   n = 0;
2006   for (i = 0; i < options->N; i++) {
2007     if (PetscCIOption(options->names[i])) continue;
2008     if (!options->used[i]) n++;
2009   }
2010   if (N) *N = n;
2011   if (names) PetscCall(PetscMalloc1(n, names));
2012   if (values) PetscCall(PetscMalloc1(n, values));
2013 
2014   n = 0;
2015   if (names || values) {
2016     for (i = 0; i < options->N; i++) {
2017       if (!options->used[i]) {
2018         if (PetscCIOption(options->names[i])) continue;
2019         if (names) (*names)[n] = options->names[i];
2020         if (values) (*values)[n] = options->values[i];
2021         n++;
2022       }
2023     }
2024   }
2025   PetscFunctionReturn(PETSC_SUCCESS);
2026 }
2027 
2028 /*@C
2029    PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using `PetscOptionsLeftGet()`.
2030 
2031    Not Collective
2032 
2033    Input Parameters:
2034 +  options - options database, use NULL for default global database
2035 .  names - names of options not used
2036 -  values - values of options not used
2037 
2038    Level: advanced
2039 
2040 .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`, `PetscOptionsLeftGet()`
2041 @*/
2042 PetscErrorCode PetscOptionsLeftRestore(PetscOptions options, PetscInt *N, char **names[], char **values[])
2043 {
2044   PetscFunctionBegin;
2045   if (N) PetscValidIntPointer(N, 2);
2046   if (names) PetscValidPointer(names, 3);
2047   if (values) PetscValidPointer(values, 4);
2048   if (N) *N = 0;
2049   if (names) PetscCall(PetscFree(*names));
2050   if (values) PetscCall(PetscFree(*values));
2051   PetscFunctionReturn(PETSC_SUCCESS);
2052 }
2053 
2054 /*@C
2055    PetscOptionsMonitorDefault - Print all options set value events using the supplied `PetscViewer`.
2056 
2057    Logically Collective
2058 
2059    Input Parameters:
2060 +  name  - option name string
2061 .  value - option value string
2062 .  source - The source for the option
2063 -  ctx - an ASCII viewer or NULL
2064 
2065    Level: intermediate
2066 
2067    Notes:
2068      If ctx is NULL, `PetscPrintf()` is used.
2069      The first MPI rank in the `PetscViewer` viewer actually prints the values, other
2070      processes may have different values set
2071 
2072      If `PetscCIEnabled` then do not print the test harness options
2073 
2074 .seealso: `PetscOptionsMonitorSet()`
2075 @*/
2076 PetscErrorCode PetscOptionsMonitorDefault(const char name[], const char value[], PetscOptionSource source, void *ctx)
2077 {
2078   PetscFunctionBegin;
2079   if (PetscCIOption(name)) PetscFunctionReturn(PETSC_SUCCESS);
2080 
2081   if (ctx) {
2082     PetscViewer viewer = (PetscViewer)ctx;
2083     if (!value) {
2084       PetscCall(PetscViewerASCIIPrintf(viewer, "Removing option: %s\n", name));
2085     } else if (!value[0]) {
2086       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
2087     } else {
2088       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
2089     }
2090   } else {
2091     MPI_Comm comm = PETSC_COMM_WORLD;
2092     if (!value) {
2093       PetscCall(PetscPrintf(comm, "Removing option: %s\n", name));
2094     } else if (!value[0]) {
2095       PetscCall(PetscPrintf(comm, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
2096     } else {
2097       PetscCall(PetscPrintf(comm, "Setting option: %s = %s (souce: %s)\n", name, value, PetscOptionSources[source]));
2098     }
2099   }
2100   PetscFunctionReturn(PETSC_SUCCESS);
2101 }
2102 
2103 /*@C
2104    PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
2105    modified the PETSc options database.
2106 
2107    Not Collective
2108 
2109    Input Parameters:
2110 +  monitor - pointer to function (if this is NULL, it turns off monitoring
2111 .  mctx    - [optional] context for private data for the
2112              monitor routine (use NULL if no context is desired)
2113 -  monitordestroy - [optional] routine that frees monitor context
2114           (may be NULL)
2115 
2116    Calling Sequence of monitor:
2117 $     monitor (const char name[], const char value[], void *mctx)
2118 
2119 +  name - option name string
2120 .  value - option value string
2121 . source - option source
2122 -  mctx  - optional monitoring context, as set by `PetscOptionsMonitorSet()`
2123 
2124    Options Database Keys:
2125    See `PetscInitialize()` for options related to option database monitoring.
2126 
2127    Notes:
2128    The default is to do nothing.  To print the name and value of options
2129    being inserted into the database, use `PetscOptionsMonitorDefault()` as the monitoring routine,
2130    with a null monitoring context.
2131 
2132    Several different monitoring routines may be set by calling
2133    `PetscOptionsMonitorSet()` multiple times; all will be called in the
2134    order in which they were set.
2135 
2136    Level: intermediate
2137 
2138 .seealso: `PetscOptionsMonitorDefault()`, `PetscInitialize()`
2139 @*/
2140 PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], PetscOptionSource, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **))
2141 {
2142   PetscOptions options = defaultoptions;
2143 
2144   PetscFunctionBegin;
2145   if (options->monitorCancel) PetscFunctionReturn(PETSC_SUCCESS);
2146   PetscCheck(options->numbermonitors < MAXOPTIONSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many PetscOptions monitors set");
2147   options->monitor[options->numbermonitors]          = monitor;
2148   options->monitordestroy[options->numbermonitors]   = monitordestroy;
2149   options->monitorcontext[options->numbermonitors++] = (void *)mctx;
2150   PetscFunctionReturn(PETSC_SUCCESS);
2151 }
2152 
2153 /*
2154    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
2155      Empty string is considered as true.
2156 */
2157 PetscErrorCode PetscOptionsStringToBool(const char value[], PetscBool *a)
2158 {
2159   PetscBool istrue, isfalse;
2160   size_t    len;
2161 
2162   PetscFunctionBegin;
2163   /* PetscStrlen() returns 0 for NULL or "" */
2164   PetscCall(PetscStrlen(value, &len));
2165   if (!len) {
2166     *a = PETSC_TRUE;
2167     PetscFunctionReturn(PETSC_SUCCESS);
2168   }
2169   PetscCall(PetscStrcasecmp(value, "TRUE", &istrue));
2170   if (istrue) {
2171     *a = PETSC_TRUE;
2172     PetscFunctionReturn(PETSC_SUCCESS);
2173   }
2174   PetscCall(PetscStrcasecmp(value, "YES", &istrue));
2175   if (istrue) {
2176     *a = PETSC_TRUE;
2177     PetscFunctionReturn(PETSC_SUCCESS);
2178   }
2179   PetscCall(PetscStrcasecmp(value, "1", &istrue));
2180   if (istrue) {
2181     *a = PETSC_TRUE;
2182     PetscFunctionReturn(PETSC_SUCCESS);
2183   }
2184   PetscCall(PetscStrcasecmp(value, "on", &istrue));
2185   if (istrue) {
2186     *a = PETSC_TRUE;
2187     PetscFunctionReturn(PETSC_SUCCESS);
2188   }
2189   PetscCall(PetscStrcasecmp(value, "FALSE", &isfalse));
2190   if (isfalse) {
2191     *a = PETSC_FALSE;
2192     PetscFunctionReturn(PETSC_SUCCESS);
2193   }
2194   PetscCall(PetscStrcasecmp(value, "NO", &isfalse));
2195   if (isfalse) {
2196     *a = PETSC_FALSE;
2197     PetscFunctionReturn(PETSC_SUCCESS);
2198   }
2199   PetscCall(PetscStrcasecmp(value, "0", &isfalse));
2200   if (isfalse) {
2201     *a = PETSC_FALSE;
2202     PetscFunctionReturn(PETSC_SUCCESS);
2203   }
2204   PetscCall(PetscStrcasecmp(value, "off", &isfalse));
2205   if (isfalse) {
2206     *a = PETSC_FALSE;
2207     PetscFunctionReturn(PETSC_SUCCESS);
2208   }
2209   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown logical value: %s", value);
2210 }
2211 
2212 /*
2213    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2214 */
2215 PetscErrorCode PetscOptionsStringToInt(const char name[], PetscInt *a)
2216 {
2217   size_t    len;
2218   PetscBool decide, tdefault, mouse;
2219 
2220   PetscFunctionBegin;
2221   PetscCall(PetscStrlen(name, &len));
2222   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
2223 
2224   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &tdefault));
2225   if (!tdefault) PetscCall(PetscStrcasecmp(name, "DEFAULT", &tdefault));
2226   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &decide));
2227   if (!decide) PetscCall(PetscStrcasecmp(name, "DECIDE", &decide));
2228   PetscCall(PetscStrcasecmp(name, "mouse", &mouse));
2229 
2230   if (tdefault) *a = PETSC_DEFAULT;
2231   else if (decide) *a = PETSC_DECIDE;
2232   else if (mouse) *a = -1;
2233   else {
2234     char *endptr;
2235     long  strtolval;
2236 
2237     strtolval = strtol(name, &endptr, 10);
2238     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);
2239 
2240 #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2241     (void)strtolval;
2242     *a = atoll(name);
2243 #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2244     (void)strtolval;
2245     *a = _atoi64(name);
2246 #else
2247     *a = (PetscInt)strtolval;
2248 #endif
2249   }
2250   PetscFunctionReturn(PETSC_SUCCESS);
2251 }
2252 
2253 #if defined(PETSC_USE_REAL___FLOAT128)
2254   #include <quadmath.h>
2255 #endif
2256 
2257 static PetscErrorCode PetscStrtod(const char name[], PetscReal *a, char **endptr)
2258 {
2259   PetscFunctionBegin;
2260 #if defined(PETSC_USE_REAL___FLOAT128)
2261   *a = strtoflt128(name, endptr);
2262 #else
2263   *a = (PetscReal)strtod(name, endptr);
2264 #endif
2265   PetscFunctionReturn(PETSC_SUCCESS);
2266 }
2267 
2268 static PetscErrorCode PetscStrtoz(const char name[], PetscScalar *a, char **endptr, PetscBool *isImaginary)
2269 {
2270   PetscBool hasi = PETSC_FALSE;
2271   char     *ptr;
2272   PetscReal strtoval;
2273 
2274   PetscFunctionBegin;
2275   PetscCall(PetscStrtod(name, &strtoval, &ptr));
2276   if (ptr == name) {
2277     strtoval = 1.;
2278     hasi     = PETSC_TRUE;
2279     if (name[0] == 'i') {
2280       ptr++;
2281     } else if (name[0] == '+' && name[1] == 'i') {
2282       ptr += 2;
2283     } else if (name[0] == '-' && name[1] == 'i') {
2284       strtoval = -1.;
2285       ptr += 2;
2286     }
2287   } else if (*ptr == 'i') {
2288     hasi = PETSC_TRUE;
2289     ptr++;
2290   }
2291   *endptr      = ptr;
2292   *isImaginary = hasi;
2293   if (hasi) {
2294 #if !defined(PETSC_USE_COMPLEX)
2295     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s contains imaginary but complex not supported ", name);
2296 #else
2297     *a = PetscCMPLX(0., strtoval);
2298 #endif
2299   } else {
2300     *a = strtoval;
2301   }
2302   PetscFunctionReturn(PETSC_SUCCESS);
2303 }
2304 
2305 /*
2306    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2307 */
2308 PetscErrorCode PetscOptionsStringToReal(const char name[], PetscReal *a)
2309 {
2310   size_t    len;
2311   PetscBool match;
2312   char     *endptr;
2313 
2314   PetscFunctionBegin;
2315   PetscCall(PetscStrlen(name, &len));
2316   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "String of length zero has no numerical value");
2317 
2318   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &match));
2319   if (!match) PetscCall(PetscStrcasecmp(name, "DEFAULT", &match));
2320   if (match) {
2321     *a = PETSC_DEFAULT;
2322     PetscFunctionReturn(PETSC_SUCCESS);
2323   }
2324 
2325   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &match));
2326   if (!match) PetscCall(PetscStrcasecmp(name, "DECIDE", &match));
2327   if (match) {
2328     *a = PETSC_DECIDE;
2329     PetscFunctionReturn(PETSC_SUCCESS);
2330   }
2331 
2332   PetscCall(PetscStrtod(name, a, &endptr));
2333   PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value", name);
2334   PetscFunctionReturn(PETSC_SUCCESS);
2335 }
2336 
2337 PetscErrorCode PetscOptionsStringToScalar(const char name[], PetscScalar *a)
2338 {
2339   PetscBool   imag1;
2340   size_t      len;
2341   PetscScalar val = 0.;
2342   char       *ptr = NULL;
2343 
2344   PetscFunctionBegin;
2345   PetscCall(PetscStrlen(name, &len));
2346   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
2347   PetscCall(PetscStrtoz(name, &val, &ptr, &imag1));
2348 #if defined(PETSC_USE_COMPLEX)
2349   if ((size_t)(ptr - name) < len) {
2350     PetscBool   imag2;
2351     PetscScalar val2;
2352 
2353     PetscCall(PetscStrtoz(ptr, &val2, &ptr, &imag2));
2354     if (imag1) PetscCheck(imag2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s: must specify imaginary component second", name);
2355     val = PetscCMPLX(PetscRealPart(val), PetscImaginaryPart(val2));
2356   }
2357 #endif
2358   PetscCheck((size_t)(ptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value ", name);
2359   *a = val;
2360   PetscFunctionReturn(PETSC_SUCCESS);
2361 }
2362 
2363 /*@C
2364    PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2365             option in the database.
2366 
2367    Not Collective
2368 
2369    Input Parameters:
2370 +  options - options database, use NULL for default global database
2371 .  pre - the string to prepend to the name or NULL
2372 -  name - the option one is seeking
2373 
2374    Output Parameters:
2375 +  ivalue - the logical value to return
2376 -  set - `PETSC_TRUE`  if found, else `PETSC_FALSE`
2377 
2378    Level: beginner
2379 
2380    Notes:
2381        TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`
2382        FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`
2383 
2384       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
2385      is equivalent to -requested_bool true
2386 
2387        If the user does not supply the option at all ivalue is NOT changed. Thus
2388      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2389 
2390 .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2391           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2392           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2393           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2394           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2395           `PetscOptionsFList()`, `PetscOptionsEList()`
2396 @*/
2397 PetscErrorCode PetscOptionsGetBool(PetscOptions options, const char pre[], const char name[], PetscBool *ivalue, PetscBool *set)
2398 {
2399   const char *value;
2400   PetscBool   flag;
2401 
2402   PetscFunctionBegin;
2403   PetscValidCharPointer(name, 3);
2404   if (ivalue) PetscValidBoolPointer(ivalue, 4);
2405   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2406   if (flag) {
2407     if (set) *set = PETSC_TRUE;
2408     PetscCall(PetscOptionsStringToBool(value, &flag));
2409     if (ivalue) *ivalue = flag;
2410   } else {
2411     if (set) *set = PETSC_FALSE;
2412   }
2413   PetscFunctionReturn(PETSC_SUCCESS);
2414 }
2415 
2416 /*@C
2417    PetscOptionsGetEList - Puts a list of option values that a single one may be selected from
2418 
2419    Not Collective
2420 
2421    Input Parameters:
2422 +  options - options database, use NULL for default global database
2423 .  pre - the string to prepend to the name or NULL
2424 .  opt - option name
2425 .  list - the possible choices (one of these must be selected, anything else is invalid)
2426 -  ntext - number of choices
2427 
2428    Output Parameters:
2429 +  value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2430 -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2431 
2432    Level: intermediate
2433 
2434    Notes:
2435     If the user does not supply the option value is NOT changed. Thus
2436      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2437 
2438    See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList`
2439 
2440 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2441           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2442           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2443           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2444           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2445           `PetscOptionsFList()`, `PetscOptionsEList()`
2446 @*/
2447 PetscErrorCode PetscOptionsGetEList(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscInt ntext, PetscInt *value, PetscBool *set)
2448 {
2449   size_t    alen, len = 0, tlen = 0;
2450   char     *svalue;
2451   PetscBool aset, flg = PETSC_FALSE;
2452   PetscInt  i;
2453 
2454   PetscFunctionBegin;
2455   PetscValidCharPointer(opt, 3);
2456   for (i = 0; i < ntext; i++) {
2457     PetscCall(PetscStrlen(list[i], &alen));
2458     if (alen > len) len = alen;
2459     tlen += len + 1;
2460   }
2461   len += 5; /* a little extra space for user mistypes */
2462   PetscCall(PetscMalloc1(len, &svalue));
2463   PetscCall(PetscOptionsGetString(options, pre, opt, svalue, len, &aset));
2464   if (aset) {
2465     PetscCall(PetscEListFind(ntext, list, svalue, value, &flg));
2466     if (!flg) {
2467       char *avail;
2468 
2469       PetscCall(PetscMalloc1(tlen, &avail));
2470       avail[0] = '\0';
2471       for (i = 0; i < ntext; i++) {
2472         PetscCall(PetscStrlcat(avail, list[i], tlen));
2473         PetscCall(PetscStrlcat(avail, " ", tlen));
2474       }
2475       PetscCall(PetscStrtolower(avail));
2476       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown option %s for -%s%s. Available options: %s", svalue, pre ? pre : "", opt + 1, avail);
2477     }
2478     if (set) *set = PETSC_TRUE;
2479   } else if (set) *set = PETSC_FALSE;
2480   PetscCall(PetscFree(svalue));
2481   PetscFunctionReturn(PETSC_SUCCESS);
2482 }
2483 
2484 /*@C
2485    PetscOptionsGetEnum - Gets the enum value for a particular option in the database.
2486 
2487    Not Collective
2488 
2489    Input Parameters:
2490 +  options - options database, use NULL for default global database
2491 .  pre - option prefix or NULL
2492 .  opt - option name
2493 -  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2494 
2495    Output Parameters:
2496 +  value - the  value to return
2497 -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2498 
2499    Level: beginner
2500 
2501    Notes:
2502     If the user does not supply the option value is NOT changed. Thus
2503      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2504 
2505           List is usually something like `PCASMTypes` or some other predefined list of enum names
2506 
2507 .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2508           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2509           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
2510           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2511           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2512           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2513           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2514 @*/
2515 PetscErrorCode PetscOptionsGetEnum(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscEnum *value, PetscBool *set)
2516 {
2517   PetscInt  ntext = 0, tval;
2518   PetscBool fset;
2519 
2520   PetscFunctionBegin;
2521   PetscValidCharPointer(opt, 3);
2522   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");
2523   PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix");
2524   ntext -= 3;
2525   PetscCall(PetscOptionsGetEList(options, pre, opt, list, ntext, &tval, &fset));
2526   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2527   if (fset) *value = (PetscEnum)tval;
2528   if (set) *set = fset;
2529   PetscFunctionReturn(PETSC_SUCCESS);
2530 }
2531 
2532 /*@C
2533    PetscOptionsGetInt - Gets the integer value for a particular option in the database.
2534 
2535    Not Collective
2536 
2537    Input Parameters:
2538 +  options - options database, use NULL for default global database
2539 .  pre - the string to prepend to the name or NULL
2540 -  name - the option one is seeking
2541 
2542    Output Parameters:
2543 +  ivalue - the integer value to return
2544 -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2545 
2546    Level: beginner
2547 
2548    Notes:
2549    If the user does not supply the option ivalue is NOT changed. Thus
2550    you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2551 
2552 .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2553           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2554           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
2555           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2556           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2557           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2558           `PetscOptionsFList()`, `PetscOptionsEList()`
2559 @*/
2560 PetscErrorCode PetscOptionsGetInt(PetscOptions options, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
2561 {
2562   const char *value;
2563   PetscBool   flag;
2564 
2565   PetscFunctionBegin;
2566   PetscValidCharPointer(name, 3);
2567   PetscValidIntPointer(ivalue, 4);
2568   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2569   if (flag) {
2570     if (!value) {
2571       if (set) *set = PETSC_FALSE;
2572     } else {
2573       if (set) *set = PETSC_TRUE;
2574       PetscCall(PetscOptionsStringToInt(value, ivalue));
2575     }
2576   } else {
2577     if (set) *set = PETSC_FALSE;
2578   }
2579   PetscFunctionReturn(PETSC_SUCCESS);
2580 }
2581 
2582 /*@C
2583    PetscOptionsGetReal - Gets the double precision value for a particular
2584    option in the database.
2585 
2586    Not Collective
2587 
2588    Input Parameters:
2589 +  options - options database, use NULL for default global database
2590 .  pre - string to prepend to each name or NULL
2591 -  name - the option one is seeking
2592 
2593    Output Parameters:
2594 +  dvalue - the double value to return
2595 -  set - `PETSC_TRUE` if found, `PETSC_FALSE` if not found
2596 
2597    Note:
2598     If the user does not supply the option dvalue is NOT changed. Thus
2599      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2600 
2601    Level: beginner
2602 
2603 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2604           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2605           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2606           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2607           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2608           `PetscOptionsFList()`, `PetscOptionsEList()`
2609 @*/
2610 PetscErrorCode PetscOptionsGetReal(PetscOptions options, const char pre[], const char name[], PetscReal *dvalue, PetscBool *set)
2611 {
2612   const char *value;
2613   PetscBool   flag;
2614 
2615   PetscFunctionBegin;
2616   PetscValidCharPointer(name, 3);
2617   PetscValidRealPointer(dvalue, 4);
2618   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2619   if (flag) {
2620     if (!value) {
2621       if (set) *set = PETSC_FALSE;
2622     } else {
2623       if (set) *set = PETSC_TRUE;
2624       PetscCall(PetscOptionsStringToReal(value, dvalue));
2625     }
2626   } else {
2627     if (set) *set = PETSC_FALSE;
2628   }
2629   PetscFunctionReturn(PETSC_SUCCESS);
2630 }
2631 
2632 /*@C
2633    PetscOptionsGetScalar - Gets the scalar value for a particular
2634    option in the database.
2635 
2636    Not Collective
2637 
2638    Input Parameters:
2639 +  options - options database, use NULL for default global database
2640 .  pre - string to prepend to each name or NULL
2641 -  name - the option one is seeking
2642 
2643    Output Parameters:
2644 +  dvalue - the double value to return
2645 -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2646 
2647    Level: beginner
2648 
2649    Usage:
2650    A complex number 2+3i must be specified with NO spaces
2651 
2652    Note:
2653     If the user does not supply the option dvalue is NOT changed. Thus
2654      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2655 
2656 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2657           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2658           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2659           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2660           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2661           `PetscOptionsFList()`, `PetscOptionsEList()`
2662 @*/
2663 PetscErrorCode PetscOptionsGetScalar(PetscOptions options, const char pre[], const char name[], PetscScalar *dvalue, PetscBool *set)
2664 {
2665   const char *value;
2666   PetscBool   flag;
2667 
2668   PetscFunctionBegin;
2669   PetscValidCharPointer(name, 3);
2670   PetscValidScalarPointer(dvalue, 4);
2671   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2672   if (flag) {
2673     if (!value) {
2674       if (set) *set = PETSC_FALSE;
2675     } else {
2676 #if !defined(PETSC_USE_COMPLEX)
2677       PetscCall(PetscOptionsStringToReal(value, dvalue));
2678 #else
2679       PetscCall(PetscOptionsStringToScalar(value, dvalue));
2680 #endif
2681       if (set) *set = PETSC_TRUE;
2682     }
2683   } else { /* flag */
2684     if (set) *set = PETSC_FALSE;
2685   }
2686   PetscFunctionReturn(PETSC_SUCCESS);
2687 }
2688 
2689 /*@C
2690    PetscOptionsGetString - Gets the string value for a particular option in
2691    the database.
2692 
2693    Not Collective
2694 
2695    Input Parameters:
2696 +  options - options database, use NULL for default global database
2697 .  pre - string to prepend to name or NULL
2698 .  name - the option one is seeking
2699 -  len - maximum length of the string including null termination
2700 
2701    Output Parameters:
2702 +  string - location to copy string
2703 -  set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2704 
2705    Level: beginner
2706 
2707    Fortran Note:
2708    The Fortran interface is slightly different from the C/C++
2709    interface (len is not used).  Sample usage in Fortran follows
2710 .vb
2711       character *20    string
2712       PetscErrorCode   ierr
2713       PetscBool        set
2714       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2715 .ve
2716 
2717    Note:
2718     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`
2719 
2720            If the user does not use the option then the string is not changed. Thus
2721            you should ALWAYS initialize the string if you access it without first checking if the set flag is true.
2722 
2723       Even if the user provided no string (for example -optionname -someotheroption) the flag is set to PETSC_TRUE (and the string is fulled with nulls).
2724 
2725 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2726           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2727           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2728           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2729           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2730           `PetscOptionsFList()`, `PetscOptionsEList()`
2731 @*/
2732 PetscErrorCode PetscOptionsGetString(PetscOptions options, const char pre[], const char name[], char string[], size_t len, PetscBool *set)
2733 {
2734   const char *value;
2735   PetscBool   flag;
2736 
2737   PetscFunctionBegin;
2738   PetscValidCharPointer(name, 3);
2739   PetscValidCharPointer(string, 4);
2740   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2741   if (!flag) {
2742     if (set) *set = PETSC_FALSE;
2743   } else {
2744     if (set) *set = PETSC_TRUE;
2745     if (value) PetscCall(PetscStrncpy(string, value, len));
2746     else PetscCall(PetscArrayzero(string, len));
2747   }
2748   PetscFunctionReturn(PETSC_SUCCESS);
2749 }
2750 
2751 char *PetscOptionsGetStringMatlab(PetscOptions options, const char pre[], const char name[])
2752 {
2753   const char *value;
2754   PetscBool   flag;
2755 
2756   PetscFunctionBegin;
2757   if (PetscOptionsFindPair(options, pre, name, &value, &flag)) PetscFunctionReturn(NULL);
2758   if (flag) PetscFunctionReturn((char *)value);
2759   PetscFunctionReturn(NULL);
2760 }
2761 
2762 /*@C
2763   PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2764   option in the database.  The values must be separated with commas with no intervening spaces.
2765 
2766   Not Collective
2767 
2768   Input Parameters:
2769 + options - options database, use NULL for default global database
2770 . pre - string to prepend to each name or NULL
2771 - name - the option one is seeking
2772 
2773   Output Parameters:
2774 + dvalue - the integer values to return
2775 . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2776 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2777 
2778   Level: beginner
2779 
2780   Note:
2781   TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE. FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
2782 
2783 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2784           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2785           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2786           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2787           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2788           `PetscOptionsFList()`, `PetscOptionsEList()`
2789 @*/
2790 PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options, const char pre[], const char name[], PetscBool dvalue[], PetscInt *nmax, PetscBool *set)
2791 {
2792   const char *svalue;
2793   char       *value;
2794   PetscInt    n = 0;
2795   PetscBool   flag;
2796   PetscToken  token;
2797 
2798   PetscFunctionBegin;
2799   PetscValidCharPointer(name, 3);
2800   PetscValidBoolPointer(dvalue, 4);
2801   PetscValidIntPointer(nmax, 5);
2802 
2803   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2804   if (!flag || !svalue) {
2805     if (set) *set = PETSC_FALSE;
2806     *nmax = 0;
2807     PetscFunctionReturn(PETSC_SUCCESS);
2808   }
2809   if (set) *set = PETSC_TRUE;
2810   PetscCall(PetscTokenCreate(svalue, ',', &token));
2811   PetscCall(PetscTokenFind(token, &value));
2812   while (value && n < *nmax) {
2813     PetscCall(PetscOptionsStringToBool(value, dvalue));
2814     PetscCall(PetscTokenFind(token, &value));
2815     dvalue++;
2816     n++;
2817   }
2818   PetscCall(PetscTokenDestroy(&token));
2819   *nmax = n;
2820   PetscFunctionReturn(PETSC_SUCCESS);
2821 }
2822 
2823 /*@C
2824   PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.
2825 
2826   Not Collective
2827 
2828   Input Parameters:
2829 + options - options database, use NULL for default global database
2830 . pre - option prefix or NULL
2831 . name - option name
2832 - list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2833 
2834   Output Parameters:
2835 + ivalue - the  enum values to return
2836 . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2837 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2838 
2839   Level: beginner
2840 
2841   Notes:
2842   The array must be passed as a comma separated list.
2843 
2844   There must be no intervening spaces between the values.
2845 
2846   list is usually something like `PCASMTypes` or some other predefined list of enum names.
2847 
2848 .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2849           `PetscOptionsGetEnum()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2850           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, `PetscOptionsName()`,
2851           `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, `PetscOptionsStringArray()`, `PetscOptionsRealArray()`,
2852           `PetscOptionsScalar()`, `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2853           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2854 @*/
2855 PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options, const char pre[], const char name[], const char *const *list, PetscEnum ivalue[], PetscInt *nmax, PetscBool *set)
2856 {
2857   const char *svalue;
2858   char       *value;
2859   PetscInt    n = 0;
2860   PetscEnum   evalue;
2861   PetscBool   flag;
2862   PetscToken  token;
2863 
2864   PetscFunctionBegin;
2865   PetscValidCharPointer(name, 3);
2866   PetscValidPointer(list, 4);
2867   PetscValidPointer(ivalue, 5);
2868   PetscValidIntPointer(nmax, 6);
2869 
2870   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2871   if (!flag || !svalue) {
2872     if (set) *set = PETSC_FALSE;
2873     *nmax = 0;
2874     PetscFunctionReturn(PETSC_SUCCESS);
2875   }
2876   if (set) *set = PETSC_TRUE;
2877   PetscCall(PetscTokenCreate(svalue, ',', &token));
2878   PetscCall(PetscTokenFind(token, &value));
2879   while (value && n < *nmax) {
2880     PetscCall(PetscEnumFind(list, value, &evalue, &flag));
2881     PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown enum value '%s' for -%s%s", svalue, pre ? pre : "", name + 1);
2882     ivalue[n++] = evalue;
2883     PetscCall(PetscTokenFind(token, &value));
2884   }
2885   PetscCall(PetscTokenDestroy(&token));
2886   *nmax = n;
2887   PetscFunctionReturn(PETSC_SUCCESS);
2888 }
2889 
2890 /*@C
2891   PetscOptionsGetIntArray - Gets an array of integer values for a particular option in the database.
2892 
2893   Not Collective
2894 
2895   Input Parameters:
2896 + options - options database, use NULL for default global database
2897 . pre - string to prepend to each name or NULL
2898 - name - the option one is seeking
2899 
2900   Output Parameters:
2901 + ivalue - the integer values to return
2902 . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
2903 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
2904 
2905   Level: beginner
2906 
2907   Notes:
2908   The array can be passed as
2909 +  a comma separated list -                                 0,1,2,3,4,5,6,7
2910 .  a range (start\-end+1) -                                 0-8
2911 .  a range with given increment (start\-end+1:inc) -        0-7:2
2912 -  a combination of values and ranges separated by commas - 0,1-8,8-15:2
2913 
2914   There must be no intervening spaces between the values.
2915 
2916 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2917           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2918           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2919           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2920           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2921           `PetscOptionsFList()`, `PetscOptionsEList()`
2922 @*/
2923 PetscErrorCode PetscOptionsGetIntArray(PetscOptions options, const char pre[], const char name[], PetscInt ivalue[], PetscInt *nmax, PetscBool *set)
2924 {
2925   const char *svalue;
2926   char       *value;
2927   PetscInt    n = 0, i, j, start, end, inc, nvalues;
2928   size_t      len;
2929   PetscBool   flag, foundrange;
2930   PetscToken  token;
2931 
2932   PetscFunctionBegin;
2933   PetscValidCharPointer(name, 3);
2934   PetscValidIntPointer(ivalue, 4);
2935   PetscValidIntPointer(nmax, 5);
2936 
2937   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2938   if (!flag || !svalue) {
2939     if (set) *set = PETSC_FALSE;
2940     *nmax = 0;
2941     PetscFunctionReturn(PETSC_SUCCESS);
2942   }
2943   if (set) *set = PETSC_TRUE;
2944   PetscCall(PetscTokenCreate(svalue, ',', &token));
2945   PetscCall(PetscTokenFind(token, &value));
2946   while (value && n < *nmax) {
2947     /* look for form  d-D where d and D are integers */
2948     foundrange = PETSC_FALSE;
2949     PetscCall(PetscStrlen(value, &len));
2950     if (value[0] == '-') i = 2;
2951     else i = 1;
2952     for (; i < (int)len; i++) {
2953       if (value[i] == '-') {
2954         PetscCheck(i != (int)len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value);
2955         value[i] = 0;
2956 
2957         PetscCall(PetscOptionsStringToInt(value, &start));
2958         inc = 1;
2959         j   = i + 1;
2960         for (; j < (int)len; j++) {
2961           if (value[j] == ':') {
2962             value[j] = 0;
2963 
2964             PetscCall(PetscOptionsStringToInt(value + j + 1, &inc));
2965             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);
2966             break;
2967           }
2968         }
2969         PetscCall(PetscOptionsStringToInt(value + i + 1, &end));
2970         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);
2971         nvalues = (end - start) / inc + (end - start) % inc;
2972         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);
2973         for (; start < end; start += inc) {
2974           *ivalue = start;
2975           ivalue++;
2976           n++;
2977         }
2978         foundrange = PETSC_TRUE;
2979         break;
2980       }
2981     }
2982     if (!foundrange) {
2983       PetscCall(PetscOptionsStringToInt(value, ivalue));
2984       ivalue++;
2985       n++;
2986     }
2987     PetscCall(PetscTokenFind(token, &value));
2988   }
2989   PetscCall(PetscTokenDestroy(&token));
2990   *nmax = n;
2991   PetscFunctionReturn(PETSC_SUCCESS);
2992 }
2993 
2994 /*@C
2995   PetscOptionsGetRealArray - Gets an array of double precision values for a
2996   particular option in the database.  The values must be separated with commas with no intervening spaces.
2997 
2998   Not Collective
2999 
3000   Input Parameters:
3001 + options - options database, use NULL for default global database
3002 . pre - string to prepend to each name or NULL
3003 - name - the option one is seeking
3004 
3005   Output Parameters:
3006 + dvalue - the double values to return
3007 . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
3008 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
3009 
3010   Level: beginner
3011 
3012 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3013           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3014           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3015           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3016           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3017           `PetscOptionsFList()`, `PetscOptionsEList()`
3018 @*/
3019 PetscErrorCode PetscOptionsGetRealArray(PetscOptions options, const char pre[], const char name[], PetscReal dvalue[], PetscInt *nmax, PetscBool *set)
3020 {
3021   const char *svalue;
3022   char       *value;
3023   PetscInt    n = 0;
3024   PetscBool   flag;
3025   PetscToken  token;
3026 
3027   PetscFunctionBegin;
3028   PetscValidCharPointer(name, 3);
3029   PetscValidRealPointer(dvalue, 4);
3030   PetscValidIntPointer(nmax, 5);
3031 
3032   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3033   if (!flag || !svalue) {
3034     if (set) *set = PETSC_FALSE;
3035     *nmax = 0;
3036     PetscFunctionReturn(PETSC_SUCCESS);
3037   }
3038   if (set) *set = PETSC_TRUE;
3039   PetscCall(PetscTokenCreate(svalue, ',', &token));
3040   PetscCall(PetscTokenFind(token, &value));
3041   while (value && n < *nmax) {
3042     PetscCall(PetscOptionsStringToReal(value, dvalue++));
3043     PetscCall(PetscTokenFind(token, &value));
3044     n++;
3045   }
3046   PetscCall(PetscTokenDestroy(&token));
3047   *nmax = n;
3048   PetscFunctionReturn(PETSC_SUCCESS);
3049 }
3050 
3051 /*@C
3052   PetscOptionsGetScalarArray - Gets an array of scalars for a
3053   particular option in the database.  The values must be separated with commas with no intervening spaces.
3054 
3055   Not Collective
3056 
3057   Input Parameters:
3058 + options - options database, use NULL for default global database
3059 . pre - string to prepend to each name or NULL
3060 - name - the option one is seeking
3061 
3062   Output Parameters:
3063 + dvalue - the scalar values to return
3064 . nmax - On input maximum number of values to retrieve, on output the actual number of values retrieved
3065 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
3066 
3067   Level: beginner
3068 
3069 .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3070           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3071           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3072           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3073           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3074           `PetscOptionsFList()`, `PetscOptionsEList()`
3075 @*/
3076 PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options, const char pre[], const char name[], PetscScalar dvalue[], PetscInt *nmax, PetscBool *set)
3077 {
3078   const char *svalue;
3079   char       *value;
3080   PetscInt    n = 0;
3081   PetscBool   flag;
3082   PetscToken  token;
3083 
3084   PetscFunctionBegin;
3085   PetscValidCharPointer(name, 3);
3086   PetscValidScalarPointer(dvalue, 4);
3087   PetscValidIntPointer(nmax, 5);
3088 
3089   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3090   if (!flag || !svalue) {
3091     if (set) *set = PETSC_FALSE;
3092     *nmax = 0;
3093     PetscFunctionReturn(PETSC_SUCCESS);
3094   }
3095   if (set) *set = PETSC_TRUE;
3096   PetscCall(PetscTokenCreate(svalue, ',', &token));
3097   PetscCall(PetscTokenFind(token, &value));
3098   while (value && n < *nmax) {
3099     PetscCall(PetscOptionsStringToScalar(value, dvalue++));
3100     PetscCall(PetscTokenFind(token, &value));
3101     n++;
3102   }
3103   PetscCall(PetscTokenDestroy(&token));
3104   *nmax = n;
3105   PetscFunctionReturn(PETSC_SUCCESS);
3106 }
3107 
3108 /*@C
3109   PetscOptionsGetStringArray - Gets an array of string values for a particular
3110   option in the database. The values must be separated with commas with no intervening spaces.
3111 
3112   Not Collective; No Fortran Support
3113 
3114   Input Parameters:
3115 + options - options database, use NULL for default global database
3116 . pre - string to prepend to name or NULL
3117 - name - the option one is seeking
3118 
3119   Output Parameters:
3120 + strings - location to copy strings
3121 . nmax - On input maximum number of strings, on output the actual number of strings found
3122 - set - `PETSC_TRUE` if found, else `PETSC_FALSE`
3123 
3124   Level: beginner
3125 
3126   Notes:
3127   The nmax parameter is used for both input and output.
3128 
3129   The user should pass in an array of pointers to char, to hold all the
3130   strings returned by this function.
3131 
3132   The user is responsible for deallocating the strings that are
3133   returned.
3134 
3135 .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
3136           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3137           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3138           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3139           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3140           `PetscOptionsFList()`, `PetscOptionsEList()`
3141 @*/
3142 PetscErrorCode PetscOptionsGetStringArray(PetscOptions options, const char pre[], const char name[], char *strings[], PetscInt *nmax, PetscBool *set)
3143 {
3144   const char *svalue;
3145   char       *value;
3146   PetscInt    n = 0;
3147   PetscBool   flag;
3148   PetscToken  token;
3149 
3150   PetscFunctionBegin;
3151   PetscValidCharPointer(name, 3);
3152   PetscValidPointer(strings, 4);
3153   PetscValidIntPointer(nmax, 5);
3154 
3155   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3156   if (!flag || !svalue) {
3157     if (set) *set = PETSC_FALSE;
3158     *nmax = 0;
3159     PetscFunctionReturn(PETSC_SUCCESS);
3160   }
3161   if (set) *set = PETSC_TRUE;
3162   PetscCall(PetscTokenCreate(svalue, ',', &token));
3163   PetscCall(PetscTokenFind(token, &value));
3164   while (value && n < *nmax) {
3165     PetscCall(PetscStrallocpy(value, &strings[n]));
3166     PetscCall(PetscTokenFind(token, &value));
3167     n++;
3168   }
3169   PetscCall(PetscTokenDestroy(&token));
3170   *nmax = n;
3171   PetscFunctionReturn(PETSC_SUCCESS);
3172 }
3173 
3174 /*@C
3175    PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one
3176 
3177    Prints a deprecation warning, unless an option is supplied to suppress.
3178 
3179    Logically Collective
3180 
3181    Input Parameters:
3182 +  pre - string to prepend to name or NULL
3183 .  oldname - the old, deprecated option
3184 .  newname - the new option, or NULL if option is purely removed
3185 .  version - a string describing the version of first deprecation, e.g. "3.9"
3186 -  info - additional information string, or NULL.
3187 
3188    Options Database Key:
3189 . -options_suppress_deprecated_warnings - do not print deprecation warnings
3190 
3191    Notes:
3192    Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
3193    Only the process of rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or
3194    `PetscObjectOptionsBegin()` prints the information
3195    If newname is provided, the old option is replaced. Otherwise, it remains
3196    in the options database.
3197    If an option is not replaced, the info argument should be used to advise the user
3198    on how to proceed.
3199    There is a limit on the length of the warning printed, so very long strings
3200    provided as info may be truncated.
3201 
3202    Level: developer
3203 
3204 .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
3205 @*/
3206 PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject, const char oldname[], const char newname[], const char version[], const char info[])
3207 {
3208   PetscBool         found, quiet;
3209   const char       *value;
3210   const char *const quietopt = "-options_suppress_deprecated_warnings";
3211   char              msg[4096];
3212   char             *prefix  = NULL;
3213   PetscOptions      options = NULL;
3214   MPI_Comm          comm    = PETSC_COMM_SELF;
3215 
3216   PetscFunctionBegin;
3217   PetscValidCharPointer(oldname, 2);
3218   PetscValidCharPointer(version, 4);
3219   if (PetscOptionsObject) {
3220     prefix  = PetscOptionsObject->prefix;
3221     options = PetscOptionsObject->options;
3222     comm    = PetscOptionsObject->comm;
3223   }
3224   PetscCall(PetscOptionsFindPair(options, prefix, oldname, &value, &found));
3225   if (found) {
3226     if (newname) {
3227       if (prefix) PetscCall(PetscOptionsPrefixPush(options, prefix));
3228       PetscCall(PetscOptionsSetValue(options, newname, value));
3229       if (prefix) PetscCall(PetscOptionsPrefixPop(options));
3230       PetscCall(PetscOptionsClearValue(options, oldname));
3231     }
3232     quiet = PETSC_FALSE;
3233     PetscCall(PetscOptionsGetBool(options, NULL, quietopt, &quiet, NULL));
3234     if (!quiet) {
3235       PetscCall(PetscStrncpy(msg, "** PETSc DEPRECATION WARNING ** : the option ", sizeof(msg)));
3236       PetscCall(PetscStrlcat(msg, oldname, sizeof(msg)));
3237       PetscCall(PetscStrlcat(msg, " is deprecated as of version ", sizeof(msg)));
3238       PetscCall(PetscStrlcat(msg, version, sizeof(msg)));
3239       PetscCall(PetscStrlcat(msg, " and will be removed in a future release.", sizeof(msg)));
3240       if (newname) {
3241         PetscCall(PetscStrlcat(msg, " Please use the option ", sizeof(msg)));
3242         PetscCall(PetscStrlcat(msg, newname, sizeof(msg)));
3243         PetscCall(PetscStrlcat(msg, " instead.", sizeof(msg)));
3244       }
3245       if (info) {
3246         PetscCall(PetscStrlcat(msg, " ", sizeof(msg)));
3247         PetscCall(PetscStrlcat(msg, info, sizeof(msg)));
3248       }
3249       PetscCall(PetscStrlcat(msg, " (Silence this warning with ", sizeof(msg)));
3250       PetscCall(PetscStrlcat(msg, quietopt, sizeof(msg)));
3251       PetscCall(PetscStrlcat(msg, ")\n", sizeof(msg)));
3252       PetscCall(PetscPrintf(comm, "%s", msg));
3253     }
3254   }
3255   PetscFunctionReturn(PETSC_SUCCESS);
3256 }
3257