Actual source code: options.c

  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() */

  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

  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: */

 12: #include <petsc/private/petscimpl.h>
 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

 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

 30: #include <petsc/private/hashtable.h>

 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: }

 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: }

 54: static inline int PetscOptEqual(const char a[], const char b[])
 55: {
 56:   return !PetscOptNameCmp(a, b);
 57: }

 59: KHASH_INIT(HO, kh_cstr_t, int, 1, PetscOptHash, PetscOptEqual)

 61: #define MAXPREFIXES        25
 62: #define MAXOPTIONSMONITORS 5

 64: const char *PetscOptionSources[] = {"code", "command line", "file", "environment"};

 66: // This table holds all the options set by the user
 67: struct _n_PetscOptions {
 68:   PetscOptions previous;

 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;

 78:   /* Hash table */
 79:   khash_t(HO) *ht;

 81:   /* Prefixes */
 82:   int  prefixind;
 83:   int  prefixstack[MAXPREFIXES];
 84:   char prefix[PETSC_MAX_OPTION_NAME];

 86:   /* Aliases */
 87:   int    Na;       /* number or aliases */
 88:   int    Naalloc;  /* number of allocated aliases */
 89:   char **aliases1; /* aliased */
 90:   char **aliases2; /* aliasee */

 92:   /* Help */
 93:   PetscBool help;       /* flag whether "-help" is in the database */
 94:   PetscBool help_intro; /* flag whether "-help intro" is in the database */

 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: };

104: static PetscOptions defaultoptions = NULL; /* the options database routines query this object for options */

106: /* list of options which precede 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: };

118: PETSC_INTERN PetscErrorCode PetscOptionsSetValue_Private(PetscOptions, const char[], const char[], int *, PetscOptionSource);
119: PETSC_INTERN PetscErrorCode PetscOptionsInsertStringYAML_Private(PetscOptions, const char[], PetscOptionSource);

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 (options->monitorFromOptions) PetscCall(PetscOptionsMonitorDefault(name, value, source, NULL));
128:   for (PetscInt i = 0; i < options->numbermonitors; i++) PetscCall((*options->monitor[i])(name, value, source, options->monitorcontext[i]));
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: /*@
133:   PetscOptionsCreate - Creates an empty options database.

135:   Logically Collective

137:   Output Parameter:
138: . options - Options database object

140:   Level: advanced

142:   Note:
143:   Though PETSc has a concept of multiple options database the current code uses a single default `PetscOptions` object

145:   Developer Notes:
146:   We may want eventually to pass a `MPI_Comm` to determine the ownership of the object

148:   This object never got developed after being introduced, it is not clear that supporting multiple `PetscOptions` objects is useful

150: .seealso: `PetscOptionsDestroy()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`
151: @*/
152: PetscErrorCode PetscOptionsCreate(PetscOptions *options)
153: {
154:   PetscFunctionBegin;
155:   PetscAssertPointer(options, 1);
156:   *options = (PetscOptions)calloc(1, sizeof(**options));
157:   PetscCheck(*options, PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate the options database");
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: /*@
162:   PetscOptionsDestroy - Destroys an option database.

164:   Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`

166:   Input Parameter:
167: . options - the `PetscOptions` object

169:   Level: advanced

171: .seealso: `PetscOptionsInsert()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsSetValue()`
172: @*/
173: PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
174: {
175:   PetscFunctionBegin;
176:   PetscAssertPointer(options, 1);
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: }

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: }

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

200:   Logically Collective

202:   Input Parameter:
203: . opt - the options obtained with `PetscOptionsCreate()`

205:   Level: advanced

207:   Notes:
208:   Use `PetscOptionsPop()` to return to the previous default options database

210:   The collectivity of this routine is complex; only the MPI ranks that call this routine will
211:   have the affect of these options. If some processes that create objects call this routine and others do
212:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
213:   on different ranks.

215:   Developer Notes:
216:   Though this functionality has been provided it has never been used in PETSc and might be removed.

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: }

229: /*@
230:   PetscOptionsPop - Pop the most recent `PetscOptionsPush()` to return to the previous default options

232:   Logically Collective on whatever communicator was associated with the call to `PetscOptionsCreate()`

234:   Level: advanced

236: .seealso: `PetscOptionsCreate()`, `PetscOptionsInsert()`, `PetscOptionsSetValue()`, `PetscOptionsLeft()`
237: @*/
238: PetscErrorCode PetscOptionsPop(void)
239: {
240:   PetscOptions current = defaultoptions;

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: }

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;

261:     PetscCall(PetscOptionsPop());
262:     PetscCall(PetscOptionsDestroy(&tmp));
263:   }
264:   PetscCall(PetscOptionsDestroy(&defaultoptions));
265:   PetscFunctionReturn(PETSC_SUCCESS);
266: }

268: /*@
269:   PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.

271:   Not Collective

273:   Input Parameter:
274: . key - string to check if valid

276:   Output Parameter:
277: . valid - `PETSC_TRUE` if a valid key

279:   Level: intermediate

281: .seealso: `PetscOptionsCreate()`, `PetscOptionsInsert()`
282: @*/
283: PetscErrorCode PetscOptionsValidKey(const char key[], PetscBool *valid)
284: {
285:   char               *ptr;
286:   PETSC_UNUSED double d;

288:   PetscFunctionBegin;
289:   if (key) PetscAssertPointer(key, 1);
290:   PetscAssertPointer(valid, 2);
291:   *valid = PETSC_FALSE;
292:   if (!key) PetscFunctionReturn(PETSC_SUCCESS);
293:   if (key[0] != '-') PetscFunctionReturn(PETSC_SUCCESS);
294:   if (key[1] == '-') key++;
295:   if (!isalpha((int)key[1])) PetscFunctionReturn(PETSC_SUCCESS);
296:   d = strtod(key, &ptr);
297:   if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) PetscFunctionReturn(PETSC_SUCCESS);
298:   *valid = PETSC_TRUE;
299:   PetscFunctionReturn(PETSC_SUCCESS);
300: }

302: static PetscErrorCode PetscOptionsInsertString_Private(PetscOptions options, const char in_str[], PetscOptionSource source)
303: {
304:   char      *first, *second;
305:   PetscToken token;

307:   PetscFunctionBegin;
308:   PetscCall(PetscTokenCreate(in_str, ' ', &token));
309:   PetscCall(PetscTokenFind(token, &first));
310:   while (first) {
311:     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;

313:     PetscCall(PetscStrcasecmp(first, "-options_file", &isfile));
314:     PetscCall(PetscStrcasecmp(first, "-options_file_yaml", &isfileyaml));
315:     PetscCall(PetscStrcasecmp(first, "-options_string_yaml", &isstringyaml));
316:     PetscCall(PetscStrcasecmp(first, "-prefix_push", &ispush));
317:     PetscCall(PetscStrcasecmp(first, "-prefix_pop", &ispop));
318:     PetscCall(PetscOptionsValidKey(first, &key));
319:     if (!key) {
320:       PetscCall(PetscTokenFind(token, &first));
321:     } else if (isfile) {
322:       PetscCall(PetscTokenFind(token, &second));
323:       PetscCall(PetscOptionsInsertFile(PETSC_COMM_SELF, options, second, PETSC_TRUE));
324:       PetscCall(PetscTokenFind(token, &first));
325:     } else if (isfileyaml) {
326:       PetscCall(PetscTokenFind(token, &second));
327:       PetscCall(PetscOptionsInsertFileYAML(PETSC_COMM_SELF, options, second, PETSC_TRUE));
328:       PetscCall(PetscTokenFind(token, &first));
329:     } else if (isstringyaml) {
330:       PetscCall(PetscTokenFind(token, &second));
331:       PetscCall(PetscOptionsInsertStringYAML_Private(options, second, source));
332:       PetscCall(PetscTokenFind(token, &first));
333:     } else if (ispush) {
334:       PetscCall(PetscTokenFind(token, &second));
335:       PetscCall(PetscOptionsPrefixPush(options, second));
336:       PetscCall(PetscTokenFind(token, &first));
337:     } else if (ispop) {
338:       PetscCall(PetscOptionsPrefixPop(options));
339:       PetscCall(PetscTokenFind(token, &first));
340:     } else {
341:       PetscCall(PetscTokenFind(token, &second));
342:       PetscCall(PetscOptionsValidKey(second, &key));
343:       if (!key) {
344:         PetscCall(PetscOptionsSetValue_Private(options, first, second, NULL, source));
345:         PetscCall(PetscTokenFind(token, &first));
346:       } else {
347:         PetscCall(PetscOptionsSetValue_Private(options, first, NULL, NULL, source));
348:         first = second;
349:       }
350:     }
351:   }
352:   PetscCall(PetscTokenDestroy(&token));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: /*@
357:   PetscOptionsInsertString - Inserts options into the database from a string

359:   Logically Collective

361:   Input Parameters:
362: + options - options object
363: - in_str  - string that contains options separated by blanks

365:   Level: intermediate

367:   The collectivity of this routine is complex; only the MPI processes that call this routine will
368:   have the affect of these options. If some processes that create objects call this routine and others do
369:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
370:   on different ranks.

372:    Contributed by Boyana Norris

374: .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
375:           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
376:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
377:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
378:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
379:           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsInsertFile()`
380: @*/
381: PetscErrorCode PetscOptionsInsertString(PetscOptions options, const char in_str[])
382: {
383:   PetscFunctionBegin;
384:   PetscCall(PetscOptionsInsertString_Private(options, in_str, PETSC_OPT_CODE));
385:   PetscFunctionReturn(PETSC_SUCCESS);
386: }

388: /*
389:     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
390: */
391: static char *Petscgetline(FILE *f)
392: {
393:   size_t size = 0;
394:   size_t len  = 0;
395:   size_t last = 0;
396:   char  *buf  = NULL;

398:   if (feof(f)) return NULL;
399:   do {
400:     size += 1024;                             /* BUFSIZ is defined as "the optimal read size for this platform" */
401:     buf = (char *)realloc((void *)buf, size); /* realloc(NULL,n) is the same as malloc(n) */
402:     /* Actually do the read. Note that fgets puts a terminal '\0' on the
403:     end of the string, so we make sure we overwrite this */
404:     if (!fgets(buf + len, 1024, f)) buf[len] = 0;
405:     PetscCallAbort(PETSC_COMM_SELF, PetscStrlen(buf, &len));
406:     last = len - 1;
407:   } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
408:   if (len) return buf;
409:   free(buf);
410:   return NULL;
411: }

413: static PetscErrorCode PetscOptionsFilename(MPI_Comm comm, const char file[], char filename[PETSC_MAX_PATH_LEN], PetscBool *yaml)
414: {
415:   char fname[PETSC_MAX_PATH_LEN + 8], path[PETSC_MAX_PATH_LEN + 8], *tail;

417:   PetscFunctionBegin;
418:   *yaml = PETSC_FALSE;
419:   PetscCall(PetscStrreplace(comm, file, fname, sizeof(fname)));
420:   PetscCall(PetscFixFilename(fname, path));
421:   PetscCall(PetscStrendswith(path, ":yaml", yaml));
422:   if (*yaml) {
423:     PetscCall(PetscStrrchr(path, ':', &tail));
424:     tail[-1] = 0; /* remove ":yaml" suffix from path */
425:   }
426:   PetscCall(PetscStrncpy(filename, path, PETSC_MAX_PATH_LEN));
427:   /* check for standard YAML and JSON filename extensions */
428:   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yaml", yaml));
429:   if (!*yaml) PetscCall(PetscStrendswith(filename, ".yml", yaml));
430:   if (!*yaml) PetscCall(PetscStrendswith(filename, ".json", yaml));
431:   if (!*yaml) { /* check file contents */
432:     PetscMPIInt rank;
433:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
434:     if (rank == 0) {
435:       FILE *fh = fopen(filename, "r");
436:       if (fh) {
437:         char buf[6] = "";
438:         if (fread(buf, 1, 6, fh) > 0) {
439:           PetscCall(PetscStrncmp(buf, "%YAML ", 6, yaml));          /* check for '%YAML' tag */
440:           if (!*yaml) PetscCall(PetscStrncmp(buf, "---", 3, yaml)); /* check for document start */
441:         }
442:         (void)fclose(fh);
443:       }
444:     }
445:     PetscCallMPI(MPI_Bcast(yaml, 1, MPIU_BOOL, 0, comm));
446:   }
447:   PetscFunctionReturn(PETSC_SUCCESS);
448: }

450: static PetscErrorCode PetscOptionsInsertFilePetsc(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
451: {
452:   char       *string, *vstring = NULL, *astring = NULL, *packed = NULL;
453:   char       *tokens[4];
454:   PetscCount  bytes;
455:   size_t      len;
456:   FILE       *fd;
457:   PetscToken  token = NULL;
458:   int         err;
459:   char       *cmatch = NULL;
460:   const char  cmt    = '#';
461:   PetscInt    line   = 1;
462:   PetscMPIInt rank, cnt = 0, acnt = 0, counts[2];
463:   PetscBool   isdir, alias = PETSC_FALSE, valid;

465:   PetscFunctionBegin;
466:   PetscCall(PetscMemzero(tokens, sizeof(tokens)));
467:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
468:   if (rank == 0) {
469:     char fpath[PETSC_MAX_PATH_LEN];
470:     char fname[PETSC_MAX_PATH_LEN];

472:     PetscCall(PetscStrreplace(PETSC_COMM_SELF, file, fpath, sizeof(fpath)));
473:     PetscCall(PetscFixFilename(fpath, fname));

475:     fd = fopen(fname, "r");
476:     PetscCall(PetscTestDirectory(fname, 'r', &isdir));
477:     PetscCheck(!isdir || !require, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified options file %s is a directory", fname);
478:     if (fd && !isdir) {
479:       PetscSegBuffer vseg, aseg;
480:       PetscCall(PetscSegBufferCreate(1, 4000, &vseg));
481:       PetscCall(PetscSegBufferCreate(1, 2000, &aseg));

483:       /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
484:       PetscCall(PetscInfo(NULL, "Opened options file %s\n", file));

486:       while ((string = Petscgetline(fd))) {
487:         /* eliminate comments from each line */
488:         PetscCall(PetscStrchr(string, cmt, &cmatch));
489:         if (cmatch) *cmatch = 0;
490:         PetscCall(PetscStrlen(string, &len));
491:         /* replace tabs, ^M, \n with " " */
492:         for (size_t i = 0; i < len; i++) {
493:           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') string[i] = ' ';
494:         }
495:         PetscCall(PetscTokenCreate(string, ' ', &token));
496:         PetscCall(PetscTokenFind(token, &tokens[0]));
497:         if (!tokens[0]) {
498:           goto destroy;
499:         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
500:           PetscCall(PetscTokenFind(token, &tokens[0]));
501:         }
502:         for (PetscInt i = 1; i < 4; i++) PetscCall(PetscTokenFind(token, &tokens[i]));
503:         if (!tokens[0]) {
504:           goto destroy;
505:         } else if (tokens[0][0] == '-') {
506:           PetscCall(PetscOptionsValidKey(tokens[0], &valid));
507:           PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": invalid option %s", fname, line, tokens[0]);
508:           PetscCall(PetscStrlen(tokens[0], &len));
509:           PetscCall(PetscSegBufferGet(vseg, len + 1, &vstring));
510:           PetscCall(PetscArraycpy(vstring, tokens[0], len));
511:           vstring[len] = ' ';
512:           if (tokens[1]) {
513:             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
514:             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]);
515:             PetscCall(PetscStrlen(tokens[1], &len));
516:             PetscCall(PetscSegBufferGet(vseg, len + 3, &vstring));
517:             vstring[0] = '"';
518:             PetscCall(PetscArraycpy(vstring + 1, tokens[1], len));
519:             vstring[len + 1] = '"';
520:             vstring[len + 2] = ' ';
521:           }
522:         } else {
523:           PetscCall(PetscStrcasecmp(tokens[0], "alias", &alias));
524:           if (alias) {
525:             PetscCall(PetscOptionsValidKey(tokens[1], &valid));
526:             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]);
527:             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]);
528:             PetscCall(PetscOptionsValidKey(tokens[2], &valid));
529:             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]);
530:             PetscCall(PetscStrlen(tokens[1], &len));
531:             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
532:             PetscCall(PetscArraycpy(astring, tokens[1], len));
533:             astring[len] = ' ';

535:             PetscCall(PetscStrlen(tokens[2], &len));
536:             PetscCall(PetscSegBufferGet(aseg, len + 1, &astring));
537:             PetscCall(PetscArraycpy(astring, tokens[2], len));
538:             astring[len] = ' ';
539:           } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown first token in options file %s line %" PetscInt_FMT ": %s", fname, line, tokens[0]);
540:         }
541:         {
542:           const char *extraToken = alias ? tokens[3] : tokens[2];
543:           PetscCheck(!extraToken, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in options file %s line %" PetscInt_FMT ": extra token %s", fname, line, extraToken);
544:         }
545:       destroy:
546:         free(string);
547:         PetscCall(PetscTokenDestroy(&token));
548:         alias = PETSC_FALSE;
549:         line++;
550:       }
551:       err = fclose(fd);
552:       PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file %s", fname);
553:       PetscCall(PetscSegBufferGetSize(aseg, &bytes)); /* size without null termination */
554:       PetscCall(PetscMPIIntCast(bytes, &acnt));
555:       PetscCall(PetscSegBufferGet(aseg, 1, &astring));
556:       astring[0] = 0;
557:       PetscCall(PetscSegBufferGetSize(vseg, &bytes)); /* size without null termination */
558:       PetscCall(PetscMPIIntCast(bytes, &cnt));
559:       PetscCall(PetscSegBufferGet(vseg, 1, &vstring));
560:       vstring[0] = 0;
561:       PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
562:       PetscCall(PetscSegBufferExtractTo(aseg, packed));
563:       PetscCall(PetscSegBufferExtractTo(vseg, packed + acnt + 1));
564:       PetscCall(PetscSegBufferDestroy(&aseg));
565:       PetscCall(PetscSegBufferDestroy(&vseg));
566:     } else PetscCheck(!require, PETSC_COMM_SELF, PETSC_ERR_USER, "Unable to open options file %s", fname);
567:   }

569:   counts[0] = acnt;
570:   counts[1] = cnt;
571:   err       = MPI_Bcast(counts, 2, MPI_INT, 0, comm);
572:   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/");
573:   acnt = counts[0];
574:   cnt  = counts[1];
575:   if (rank) PetscCall(PetscMalloc1(2 + acnt + cnt, &packed));
576:   if (acnt || cnt) {
577:     PetscCallMPI(MPI_Bcast(packed, 2 + acnt + cnt, MPI_CHAR, 0, comm));
578:     astring = packed;
579:     vstring = packed + acnt + 1;
580:   }

582:   if (acnt) {
583:     PetscCall(PetscTokenCreate(astring, ' ', &token));
584:     PetscCall(PetscTokenFind(token, &tokens[0]));
585:     while (tokens[0]) {
586:       PetscCall(PetscTokenFind(token, &tokens[1]));
587:       PetscCall(PetscOptionsSetAlias(options, tokens[0], tokens[1]));
588:       PetscCall(PetscTokenFind(token, &tokens[0]));
589:     }
590:     PetscCall(PetscTokenDestroy(&token));
591:   }

593:   if (cnt) PetscCall(PetscOptionsInsertString_Private(options, vstring, PETSC_OPT_FILE));
594:   PetscCall(PetscFree(packed));
595:   PetscFunctionReturn(PETSC_SUCCESS);
596: }

598: /*@
599:   PetscOptionsInsertFile - Inserts options into the database from a file.

601:   Collective

603:   Input Parameters:
604: + comm    - the processes that will share the options (usually `PETSC_COMM_WORLD`)
605: . options - options database, use `NULL` for default global database
606: . file    - name of file,
607:            ".yml" and ".yaml" filename extensions are inserted as YAML options,
608:            append ":yaml" to filename to force YAML options.
609: - require - if `PETSC_TRUE` will generate an error if the file does not exist

611:   Level: developer

613:   Notes:
614:   Use  # for lines that are comments and which should be ignored.
615:   Usually, instead of using this command, one should list the file name in the call to `PetscInitialize()`, this insures that certain options
616:   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
617:   calls to `XXXSetFromOptions()`, it should not be used for options listed under PetscInitialize().
618:   The collectivity of this routine is complex; only the MPI processes in comm will
619:   have the effect of these options. If some processes that create objects call this routine and others do
620:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
621:   on different ranks.

623: .seealso: `PetscOptionsSetValue()`, `PetscOptionsView()`, `PetscOptionsHasName()`, `PetscOptionsGetInt()`,
624:           `PetscOptionsGetReal()`, `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
625:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
626:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
627:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
628:           `PetscOptionsFList()`, `PetscOptionsEList()`
629: @*/
630: PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm, PetscOptions options, const char file[], PetscBool require)
631: {
632:   char      filename[PETSC_MAX_PATH_LEN];
633:   PetscBool yaml;

635:   PetscFunctionBegin;
636:   PetscCall(PetscOptionsFilename(comm, file, filename, &yaml));
637:   if (yaml) {
638:     PetscCall(PetscOptionsInsertFileYAML(comm, options, filename, require));
639:   } else {
640:     PetscCall(PetscOptionsInsertFilePetsc(comm, options, filename, require));
641:   }
642:   PetscFunctionReturn(PETSC_SUCCESS);
643: }

645: /*@C
646:   PetscOptionsInsertArgs - Inserts options into the database from a array of strings

648:   Logically Collective

650:   Input Parameters:
651: + options - options object
652: . argc    - the array length
653: - args    - the string array

655:   Level: intermediate

657: .seealso: `PetscOptions`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`
658: @*/
659: PetscErrorCode PetscOptionsInsertArgs(PetscOptions options, int argc, char *args[])
660: {
661:   MPI_Comm     comm  = PETSC_COMM_WORLD;
662:   int          left  = PetscMax(argc, 0);
663:   char *const *eargs = args;

665:   PetscFunctionBegin;
666:   while (left) {
667:     PetscBool isfile, isfileyaml, isstringyaml, ispush, ispop, key;
668:     PetscCall(PetscStrcasecmp(eargs[0], "-options_file", &isfile));
669:     PetscCall(PetscStrcasecmp(eargs[0], "-options_file_yaml", &isfileyaml));
670:     PetscCall(PetscStrcasecmp(eargs[0], "-options_string_yaml", &isstringyaml));
671:     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_push", &ispush));
672:     PetscCall(PetscStrcasecmp(eargs[0], "-prefix_pop", &ispop));
673:     PetscCall(PetscOptionsValidKey(eargs[0], &key));
674:     if (!key) {
675:       eargs++;
676:       left--;
677:     } else if (isfile) {
678:       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file filename option");
679:       PetscCall(PetscOptionsInsertFile(comm, options, eargs[1], PETSC_TRUE));
680:       eargs += 2;
681:       left -= 2;
682:     } else if (isfileyaml) {
683:       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing filename for -options_file_yaml filename option");
684:       PetscCall(PetscOptionsInsertFileYAML(comm, options, eargs[1], PETSC_TRUE));
685:       eargs += 2;
686:       left -= 2;
687:     } else if (isstringyaml) {
688:       PetscCheck(left > 1 && eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing string for -options_string_yaml string option");
689:       PetscCall(PetscOptionsInsertStringYAML_Private(options, eargs[1], PETSC_OPT_CODE));
690:       eargs += 2;
691:       left -= 2;
692:     } else if (ispush) {
693:       PetscCheck(left > 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option");
694:       PetscCheck(eargs[1][0] != '-', PETSC_COMM_SELF, PETSC_ERR_USER, "Missing prefix for -prefix_push option (prefixes cannot start with '-')");
695:       PetscCall(PetscOptionsPrefixPush(options, eargs[1]));
696:       eargs += 2;
697:       left -= 2;
698:     } else if (ispop) {
699:       PetscCall(PetscOptionsPrefixPop(options));
700:       eargs++;
701:       left--;
702:     } else {
703:       PetscBool nextiskey = PETSC_FALSE;
704:       if (left >= 2) PetscCall(PetscOptionsValidKey(eargs[1], &nextiskey));
705:       if (left < 2 || nextiskey) {
706:         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], NULL, NULL, PETSC_OPT_COMMAND_LINE));
707:         eargs++;
708:         left--;
709:       } else {
710:         PetscCall(PetscOptionsSetValue_Private(options, eargs[0], eargs[1], NULL, PETSC_OPT_COMMAND_LINE));
711:         eargs += 2;
712:         left -= 2;
713:       }
714:     }
715:   }
716:   PetscFunctionReturn(PETSC_SUCCESS);
717: }

719: static inline PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt, const char *val[], const PetscBool set[], PetscBool *flg)
720: {
721:   PetscFunctionBegin;
722:   if (set[opt]) {
723:     PetscCall(PetscOptionsStringToBool(val[opt], flg));
724:   } else *flg = PETSC_FALSE;
725:   PetscFunctionReturn(PETSC_SUCCESS);
726: }

728: /* Process options with absolute precedence, these are only processed from the command line, not the environment or files */
729: static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options, int argc, char *args[], PetscBool *skip_petscrc, PetscBool *skip_petscrc_set)
730: {
731:   const char *const *opt = precedentOptions;
732:   const size_t       n   = PO_NUM;
733:   size_t             o;
734:   int                a;
735:   const char       **val;
736:   char             **cval;
737:   PetscBool         *set, unneeded;

739:   PetscFunctionBegin;
740:   PetscCall(PetscCalloc2(n, &cval, n, &set));
741:   val = (const char **)cval;

743:   /* Look for options possibly set using PetscOptionsSetValue beforehand */
744:   for (o = 0; o < n; o++) PetscCall(PetscOptionsFindPair(options, NULL, opt[o], &val[o], &set[o]));

746:   /* Loop through all args to collect last occurring value of each option */
747:   for (a = 1; a < argc; a++) {
748:     PetscBool valid, eq;

750:     PetscCall(PetscOptionsValidKey(args[a], &valid));
751:     if (!valid) continue;
752:     for (o = 0; o < n; o++) {
753:       PetscCall(PetscStrcasecmp(args[a], opt[o], &eq));
754:       if (eq) {
755:         set[o] = PETSC_TRUE;
756:         if (a == argc - 1 || !args[a + 1] || !args[a + 1][0] || args[a + 1][0] == '-') val[o] = NULL;
757:         else val[o] = args[a + 1];
758:         break;
759:       }
760:     }
761:   }

763:   /* Process flags */
764:   PetscCall(PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro));
765:   if (options->help_intro) options->help = PETSC_TRUE;
766:   else PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_HELP, val, set, &options->help));
767:   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_CI_ENABLE, val, set, &unneeded));
768:   /* need to manage PO_CI_ENABLE option before the PetscOptionsMonitor is turned on, so its setting is not monitored */
769:   if (set[PO_CI_ENABLE]) PetscCall(PetscOptionsSetValue_Private(options, opt[PO_CI_ENABLE], val[PO_CI_ENABLE], &a, PETSC_OPT_COMMAND_LINE));
770:   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL, val, set, &options->monitorCancel));
771:   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR, val, set, &options->monitorFromOptions));
772:   PetscCall(PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC, val, set, skip_petscrc));
773:   *skip_petscrc_set = set[PO_SKIP_PETSCRC];

775:   /* Store precedent options in database and mark them as used */
776:   for (o = 1; o < n; o++) {
777:     if (set[o]) {
778:       PetscCall(PetscOptionsSetValue_Private(options, opt[o], val[o], &a, PETSC_OPT_COMMAND_LINE));
779:       options->used[a] = PETSC_TRUE;
780:     }
781:   }
782:   PetscCall(PetscFree2(cval, set));
783:   options->precedentProcessed = PETSC_TRUE;
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

787: static inline PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options, const char name[], PetscBool *flg)
788: {
789:   PetscFunctionBegin;
790:   PetscAssertPointer(flg, 3);
791:   *flg = PETSC_FALSE;
792:   if (options->precedentProcessed) {
793:     for (int i = 0; i < PO_NUM; ++i) {
794:       if (!PetscOptNameCmp(precedentOptions[i], name)) {
795:         /* check if precedent option has been set already */
796:         PetscCall(PetscOptionsFindPair(options, NULL, name, NULL, flg));
797:         if (*flg) break;
798:       }
799:     }
800:   }
801:   PetscFunctionReturn(PETSC_SUCCESS);
802: }

804: /*@C
805:   PetscOptionsInsert - Inserts into the options database from the command line,
806:   the environmental variable and a file.

808:   Collective on `PETSC_COMM_WORLD`

810:   Input Parameters:
811: + options - options database or `NULL` for the default global database
812: . argc    - count of number of command line arguments
813: . args    - the command line arguments
814: - file    - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
815:           Use `NULL` or empty string to not check for code specific file.
816:           Also checks ~/.petscrc, .petscrc and petscrc.
817:           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.

819:   Options Database Keys:
820: + -options_file <filename>      - read options from a file
821: - -options_file_yaml <filename> - read options from a YAML file

823:   Level: advanced

825:   Notes:
826:   Since `PetscOptionsInsert()` is automatically called by `PetscInitialize()`,
827:   the user does not typically need to call this routine. `PetscOptionsInsert()`
828:   can be called several times, adding additional entries into the database.

830:   See `PetscInitialize()` for options related to option database monitoring.

832: .seealso: `PetscOptionsDestroy()`, `PetscOptionsView()`, `PetscOptionsInsertString()`, `PetscOptionsInsertFile()`,
833:           `PetscInitialize()`
834: @*/
835: PetscErrorCode PetscOptionsInsert(PetscOptions options, int *argc, char ***args, const char file[])
836: {
837:   MPI_Comm    comm = PETSC_COMM_WORLD;
838:   PetscMPIInt rank;
839:   PetscBool   hasArgs     = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
840:   PetscBool   skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;
841:   char       *eoptions = NULL;
842:   size_t      len      = 0;

844:   PetscFunctionBegin;
845:   PetscCheck(!hasArgs || (args && *args), comm, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
846:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

848:   if (!options) {
849:     PetscCall(PetscOptionsCreateDefault());
850:     options = defaultoptions;
851:   }
852:   if (hasArgs) {
853:     /* process options with absolute precedence */
854:     PetscCall(PetscOptionsProcessPrecedentFlags(options, *argc, *args, &skipPetscrc, &skipPetscrcSet));
855:     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci", &PetscCIEnabled, NULL));
856:   }
857:   if (file && file[0]) {
858:     PetscCall(PetscOptionsInsertFile(comm, options, file, PETSC_TRUE));
859:     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
860:     if (!skipPetscrcSet) PetscCall(PetscOptionsGetBool(options, NULL, "-skip_petscrc", &skipPetscrc, NULL));
861:   }
862:   if (!skipPetscrc) {
863:     char filename[PETSC_MAX_PATH_LEN];

865:     PetscCall(PetscGetHomeDirectory(filename, sizeof(filename)));
866:     PetscCallMPI(MPI_Bcast(filename, (int)sizeof(filename), MPI_CHAR, 0, comm));
867:     if (filename[0]) PetscCall(PetscStrlcat(filename, "/.petscrc", sizeof(filename)));
868:     PetscCall(PetscOptionsInsertFile(comm, options, filename, PETSC_FALSE));
869:     PetscCall(PetscOptionsInsertFile(comm, options, ".petscrc", PETSC_FALSE));
870:     PetscCall(PetscOptionsInsertFile(comm, options, "petscrc", PETSC_FALSE));
871:   }

873:   /* insert environment options */
874:   if (rank == 0) {
875:     eoptions = (char *)getenv("PETSC_OPTIONS");
876:     PetscCall(PetscStrlen(eoptions, &len));
877:   }
878:   PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
879:   if (len) {
880:     if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
881:     PetscCallMPI(MPI_Bcast(eoptions, (PetscMPIInt)len, MPI_CHAR, 0, comm));
882:     if (rank) eoptions[len] = 0;
883:     PetscCall(PetscOptionsInsertString_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
884:     if (rank) PetscCall(PetscFree(eoptions));
885:   }

887:   /* insert YAML environment options */
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, (PetscMPIInt)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:   }

901:   /* insert command line options here because they take precedence over arguments in petscrc/environment */
902:   if (hasArgs) PetscCall(PetscOptionsInsertArgs(options, *argc - 1, *args + 1));
903:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-petsc_ci_portable_error_output", &PetscCIEnabledPortableErrorOutput, NULL));
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }

907: /* These options are not printed with PetscOptionsView() or PetscOptionsMonitor() when PetscCIEnabled is on */
908: /* TODO: get the list from the test harness, do not have it hardwired here. Maybe from gmakegentest.py */
909: static const char *PetscCIOptions[] = {"malloc_debug", "malloc_dump", "malloc_test", "malloc", "nox", "nox_warning", "display", "saws_port_auto_select", "saws_port_auto_select_silent", "vecscatter_mpi1", "check_pointer_intensity", "cuda_initialize", "error_output_stdout", "use_gpu_aware_mpi", "checkfunctionlist", "fp_trap", "petsc_ci", "petsc_ci_portable_error_output", "options_left"};

911: static PetscBool PetscCIOption(const char *name)
912: {
913:   PetscInt  idx;
914:   PetscBool found;

916:   if (!PetscCIEnabled) return PETSC_FALSE;
917:   PetscCallAbort(PETSC_COMM_SELF, PetscEListFind(PETSC_STATIC_ARRAY_LENGTH(PetscCIOptions), PetscCIOptions, name, &idx, &found));
918:   return found;
919: }

921: /*@
922:   PetscOptionsView - Prints the options that have been loaded. This is
923:   useful for debugging purposes.

925:   Logically Collective, No Fortran Support

927:   Input Parameters:
928: + options - options database, use `NULL` for default global database
929: - viewer  - must be an `PETSCVIEWERASCII` viewer

931:   Options Database Key:
932: . -options_view - Activates `PetscOptionsView()` within `PetscFinalize()`

934:   Level: advanced

936:   Note:
937:   Only the MPI rank 0 of the `MPI_Comm` used to create view prints the option values. Other processes
938:   may have different values but they are not printed.

940: .seealso: `PetscOptionsAllUsed()`
941: @*/
942: PetscErrorCode PetscOptionsView(PetscOptions options, PetscViewer viewer)
943: {
944:   PetscInt  i, N = 0;
945:   PetscBool isascii;

947:   PetscFunctionBegin;
949:   options = options ? options : defaultoptions;
950:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
951:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
952:   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Only supports ASCII viewer");

954:   for (i = 0; i < options->N; i++) {
955:     if (PetscCIOption(options->names[i])) continue;
956:     N++;
957:   }

959:   if (!N) {
960:     PetscCall(PetscViewerASCIIPrintf(viewer, "#No PETSc Option Table entries\n"));
961:     PetscFunctionReturn(PETSC_SUCCESS);
962:   }

964:   PetscCall(PetscViewerASCIIPrintf(viewer, "#PETSc Option Table entries:\n"));
965:   for (i = 0; i < options->N; i++) {
966:     if (PetscCIOption(options->names[i])) continue;
967:     if (options->values[i]) {
968:       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s %s", options->names[i], options->values[i]));
969:     } else {
970:       PetscCall(PetscViewerASCIIPrintf(viewer, "-%s", options->names[i]));
971:     }
972:     PetscCall(PetscViewerASCIIPrintf(viewer, " # (source: %s)\n", PetscOptionSources[options->source[i]]));
973:   }
974:   PetscCall(PetscViewerASCIIPrintf(viewer, "#End of PETSc Option Table entries\n"));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }

978: /*
979:    Called by error handlers to print options used in run
980: */
981: PetscErrorCode PetscOptionsLeftError(void)
982: {
983:   PetscInt i, nopt = 0;

985:   for (i = 0; i < defaultoptions->N; i++) {
986:     if (!defaultoptions->used[i]) {
987:       if (PetscCIOption(defaultoptions->names[i])) continue;
988:       nopt++;
989:     }
990:   }
991:   if (nopt) {
992:     PetscCall((*PetscErrorPrintf)("WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!\n"));
993:     for (i = 0; i < defaultoptions->N; i++) {
994:       if (!defaultoptions->used[i]) {
995:         if (PetscCIOption(defaultoptions->names[i])) continue;
996:         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]]));
997:         else PetscCall((*PetscErrorPrintf)("  Option left: name:-%s (no value) source: %s\n", defaultoptions->names[i], PetscOptionSources[defaultoptions->source[i]]));
998:       }
999:     }
1000:   }
1001:   return PETSC_SUCCESS;
1002: }

1004: PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
1005: {
1006:   PetscInt     i, N = 0;
1007:   PetscOptions options = defaultoptions;

1009:   for (i = 0; i < options->N; i++) {
1010:     if (PetscCIOption(options->names[i])) continue;
1011:     N++;
1012:   }

1014:   if (N) {
1015:     PetscCall((*PetscErrorPrintf)("PETSc Option Table entries:\n"));
1016:   } else {
1017:     PetscCall((*PetscErrorPrintf)("No PETSc Option Table entries\n"));
1018:   }
1019:   for (i = 0; i < options->N; i++) {
1020:     if (PetscCIOption(options->names[i])) continue;
1021:     if (options->values[i]) {
1022:       PetscCall((*PetscErrorPrintf)("-%s %s (source: %s)\n", options->names[i], options->values[i], PetscOptionSources[options->source[i]]));
1023:     } else {
1024:       PetscCall((*PetscErrorPrintf)("-%s (source: %s)\n", options->names[i], PetscOptionSources[options->source[i]]));
1025:     }
1026:   }
1027:   return PETSC_SUCCESS;
1028: }

1030: /*@
1031:   PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.

1033:   Logically Collective

1035:   Input Parameters:
1036: + options - options database, or `NULL` for the default global database
1037: - prefix  - The string to append to the existing prefix

1039:   Options Database Keys:
1040: + -prefix_push <some_prefix_> - push the given prefix
1041: - -prefix_pop                 - pop the last prefix

1043:   Level: advanced

1045:   Notes:
1046:   It is common to use this in conjunction with `-options_file` as in
1047: .vb
1048:  -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
1049: .ve
1050:   where the files no longer require all options to be prefixed with `-system2_`.

1052:   The collectivity of this routine is complex; only the MPI processes that call this routine will
1053:   have the affect of these options. If some processes that create objects call this routine and others do
1054:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1055:   on different ranks.

1057: .seealso: `PetscOptionsPrefixPop()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
1058: @*/
1059: PetscErrorCode PetscOptionsPrefixPush(PetscOptions options, const char prefix[])
1060: {
1061:   size_t    n;
1062:   PetscInt  start;
1063:   char      key[PETSC_MAX_OPTION_NAME + 1];
1064:   PetscBool valid;

1066:   PetscFunctionBegin;
1067:   PetscAssertPointer(prefix, 2);
1068:   options = options ? options : defaultoptions;
1069:   PetscCheck(options->prefixind < MAXPREFIXES, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum depth of prefix stack %d exceeded, recompile src/sys/objects/options.c with larger value for MAXPREFIXES", MAXPREFIXES);
1070:   key[0] = '-'; /* keys must start with '-' */
1071:   PetscCall(PetscStrncpy(key + 1, prefix, sizeof(key) - 1));
1072:   PetscCall(PetscOptionsValidKey(key, &valid));
1073:   if (!valid && options->prefixind > 0 && isdigit((int)prefix[0])) valid = PETSC_TRUE; /* If the prefix stack is not empty, make numbers a valid prefix */
1074:   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" : "");
1075:   start = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
1076:   PetscCall(PetscStrlen(prefix, &n));
1077:   PetscCheck(n + 1 <= sizeof(options->prefix) - start, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Maximum prefix length %zu exceeded", sizeof(options->prefix));
1078:   PetscCall(PetscArraycpy(options->prefix + start, prefix, n + 1));
1079:   options->prefixstack[options->prefixind++] = (int)(start + n);
1080:   PetscFunctionReturn(PETSC_SUCCESS);
1081: }

1083: /*@
1084:   PetscOptionsPrefixPop - Remove the latest options prefix, see `PetscOptionsPrefixPush()` for details

1086:   Logically Collective on the `MPI_Comm` used when called `PetscOptionsPrefixPush()`

1088:   Input Parameter:
1089: . options - options database, or `NULL` for the default global database

1091:   Level: advanced

1093: .seealso: `PetscOptionsPrefixPush()`, `PetscOptionsPush()`, `PetscOptionsPop()`, `PetscOptionsCreate()`, `PetscOptionsSetValue()`
1094: @*/
1095: PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1096: {
1097:   PetscInt offset;

1099:   PetscFunctionBegin;
1100:   options = options ? options : defaultoptions;
1101:   PetscCheck(options->prefixind >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "More prefixes popped than pushed");
1102:   options->prefixind--;
1103:   offset                  = options->prefixind ? options->prefixstack[options->prefixind - 1] : 0;
1104:   options->prefix[offset] = 0;
1105:   PetscFunctionReturn(PETSC_SUCCESS);
1106: }

1108: /*@
1109:   PetscOptionsClear - Removes all options form the database leaving it empty.

1111:   Logically Collective

1113:   Input Parameter:
1114: . options - options database, use `NULL` for the default global database

1116:   Level: developer

1118:   Note:
1119:   The collectivity of this routine is complex; only the MPI processes that call this routine will
1120:   have the affect of these options. If some processes that create objects call this routine and others do
1121:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1122:   on different ranks.

1124:   Developer Note:
1125:   Uses `free()` directly because the current option values were set with `malloc()`

1127: .seealso: `PetscOptionsInsert()`
1128: @*/
1129: PetscErrorCode PetscOptionsClear(PetscOptions options)
1130: {
1131:   PetscInt i;

1133:   PetscFunctionBegin;
1134:   options = options ? options : defaultoptions;
1135:   if (!options) PetscFunctionReturn(PETSC_SUCCESS);

1137:   for (i = 0; i < options->N; i++) {
1138:     if (options->names[i]) free(options->names[i]);
1139:     if (options->values[i]) free(options->values[i]);
1140:   }
1141:   options->N = 0;
1142:   free(options->names);
1143:   free(options->values);
1144:   free(options->used);
1145:   free(options->source);
1146:   options->names  = NULL;
1147:   options->values = NULL;
1148:   options->used   = NULL;
1149:   options->source = NULL;
1150:   options->Nalloc = 0;

1152:   for (i = 0; i < options->Na; i++) {
1153:     free(options->aliases1[i]);
1154:     free(options->aliases2[i]);
1155:   }
1156:   options->Na = 0;
1157:   free(options->aliases1);
1158:   free(options->aliases2);
1159:   options->aliases1 = options->aliases2 = NULL;
1160:   options->Naalloc                      = 0;

1162:   /* destroy hash table */
1163:   kh_destroy(HO, options->ht);
1164:   options->ht = NULL;

1166:   options->prefixind  = 0;
1167:   options->prefix[0]  = 0;
1168:   options->help       = PETSC_FALSE;
1169:   options->help_intro = PETSC_FALSE;
1170:   PetscFunctionReturn(PETSC_SUCCESS);
1171: }

1173: /*@
1174:   PetscOptionsSetAlias - Makes a key and alias for another key

1176:   Logically Collective

1178:   Input Parameters:
1179: + options - options database, or `NULL` for default global database
1180: . newname - the alias
1181: - oldname - the name that alias will refer to

1183:   Level: advanced

1185:   Note:
1186:   The collectivity of this routine is complex; only the MPI processes that call this routine will
1187:   have the affect of these options. If some processes that create objects call this routine and others do
1188:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1189:   on different ranks.

1191:   Developer Note:
1192:   Uses `malloc()` directly because PETSc may not be initialized yet.

1194: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`,
1195:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1196:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1197:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1198:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1199:           `PetscOptionsFList()`, `PetscOptionsEList()`
1200: @*/
1201: PetscErrorCode PetscOptionsSetAlias(PetscOptions options, const char newname[], const char oldname[])
1202: {
1203:   size_t    len;
1204:   PetscBool valid;

1206:   PetscFunctionBegin;
1207:   PetscAssertPointer(newname, 2);
1208:   PetscAssertPointer(oldname, 3);
1209:   options = options ? options : defaultoptions;
1210:   PetscCall(PetscOptionsValidKey(newname, &valid));
1211:   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliased option %s", newname);
1212:   PetscCall(PetscOptionsValidKey(oldname, &valid));
1213:   PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid aliasee option %s", oldname);

1215:   if (options->Na == options->Naalloc) {
1216:     char **tmpA1, **tmpA2;

1218:     options->Naalloc = PetscMax(4, options->Naalloc * 2);
1219:     tmpA1            = (char **)malloc(options->Naalloc * sizeof(char *));
1220:     tmpA2            = (char **)malloc(options->Naalloc * sizeof(char *));
1221:     for (int i = 0; i < options->Na; ++i) {
1222:       tmpA1[i] = options->aliases1[i];
1223:       tmpA2[i] = options->aliases2[i];
1224:     }
1225:     free(options->aliases1);
1226:     free(options->aliases2);
1227:     options->aliases1 = tmpA1;
1228:     options->aliases2 = tmpA2;
1229:   }
1230:   newname++;
1231:   oldname++;
1232:   PetscCall(PetscStrlen(newname, &len));
1233:   options->aliases1[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1234:   PetscCall(PetscStrncpy(options->aliases1[options->Na], newname, len + 1));
1235:   PetscCall(PetscStrlen(oldname, &len));
1236:   options->aliases2[options->Na] = (char *)malloc((len + 1) * sizeof(char));
1237:   PetscCall(PetscStrncpy(options->aliases2[options->Na], oldname, len + 1));
1238:   ++options->Na;
1239:   PetscFunctionReturn(PETSC_SUCCESS);
1240: }

1242: /*@
1243:   PetscOptionsSetValue - Sets an option name-value pair in the options
1244:   database, overriding whatever is already present.

1246:   Logically Collective

1248:   Input Parameters:
1249: + options - options database, use `NULL` for the default global database
1250: . name    - name of option, this SHOULD have the - prepended
1251: - value   - the option value (not used for all options, so can be `NULL`)

1253:   Level: intermediate

1255:   Note:
1256:   This function can be called BEFORE `PetscInitialize()`

1258:   The collectivity of this routine is complex; only the MPI processes that call this routine will
1259:   have the affect of these options. If some processes that create objects call this routine and others do
1260:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1261:   on different ranks.

1263:   Developer Note:
1264:   Uses `malloc()` directly because PETSc may not be initialized yet.

1266: .seealso: `PetscOptionsInsert()`, `PetscOptionsClearValue()`
1267: @*/
1268: PetscErrorCode PetscOptionsSetValue(PetscOptions options, const char name[], const char value[])
1269: {
1270:   PetscFunctionBegin;
1271:   PetscCall(PetscOptionsSetValue_Private(options, name, value, NULL, PETSC_OPT_CODE));
1272:   PetscFunctionReturn(PETSC_SUCCESS);
1273: }

1275: PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options, const char name[], const char value[], int *pos, PetscOptionSource source)
1276: {
1277:   size_t    len;
1278:   int       n, i;
1279:   char    **names;
1280:   char      fullname[PETSC_MAX_OPTION_NAME] = "";
1281:   PetscBool flg;

1283:   PetscFunctionBegin;
1284:   if (!options) {
1285:     PetscCall(PetscOptionsCreateDefault());
1286:     options = defaultoptions;
1287:   }
1288:   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "name %s must start with '-'", name);

1290:   PetscCall(PetscOptionsSkipPrecedent(options, name, &flg));
1291:   if (flg) PetscFunctionReturn(PETSC_SUCCESS);

1293:   name++; /* skip starting dash */

1295:   if (options->prefixind > 0) {
1296:     strncpy(fullname, options->prefix, sizeof(fullname));
1297:     fullname[sizeof(fullname) - 1] = 0;
1298:     strncat(fullname, name, sizeof(fullname) - strlen(fullname) - 1);
1299:     fullname[sizeof(fullname) - 1] = 0;
1300:     name                           = fullname;
1301:   }

1303:   /* check against aliases */
1304:   for (i = 0; i < options->Na; i++) {
1305:     int result = PetscOptNameCmp(options->aliases1[i], name);
1306:     if (!result) {
1307:       name = options->aliases2[i];
1308:       break;
1309:     }
1310:   }

1312:   /* slow search */
1313:   n     = options->N;
1314:   names = options->names;
1315:   for (i = 0; i < options->N; i++) {
1316:     int result = PetscOptNameCmp(names[i], name);
1317:     if (!result) {
1318:       n = i;
1319:       goto setvalue;
1320:     } else if (result > 0) {
1321:       n = i;
1322:       break;
1323:     }
1324:   }
1325:   if (options->N == options->Nalloc) {
1326:     char             **names, **values;
1327:     PetscBool         *used;
1328:     PetscOptionSource *source;

1330:     options->Nalloc = PetscMax(10, options->Nalloc * 2);
1331:     names           = (char **)malloc(options->Nalloc * sizeof(char *));
1332:     values          = (char **)malloc(options->Nalloc * sizeof(char *));
1333:     used            = (PetscBool *)malloc(options->Nalloc * sizeof(PetscBool));
1334:     source          = (PetscOptionSource *)malloc(options->Nalloc * sizeof(PetscOptionSource));
1335:     for (int i = 0; i < options->N; ++i) {
1336:       names[i]  = options->names[i];
1337:       values[i] = options->values[i];
1338:       used[i]   = options->used[i];
1339:       source[i] = options->source[i];
1340:     }
1341:     free(options->names);
1342:     free(options->values);
1343:     free(options->used);
1344:     free(options->source);
1345:     options->names  = names;
1346:     options->values = values;
1347:     options->used   = used;
1348:     options->source = source;
1349:   }

1351:   /* shift remaining values up 1 */
1352:   for (i = options->N; i > n; i--) {
1353:     options->names[i]  = options->names[i - 1];
1354:     options->values[i] = options->values[i - 1];
1355:     options->used[i]   = options->used[i - 1];
1356:     options->source[i] = options->source[i - 1];
1357:   }
1358:   options->names[n]  = NULL;
1359:   options->values[n] = NULL;
1360:   options->used[n]   = PETSC_FALSE;
1361:   options->source[n] = PETSC_OPT_CODE;
1362:   options->N++;

1364:   /* destroy hash table */
1365:   kh_destroy(HO, options->ht);
1366:   options->ht = NULL;

1368:   /* set new name */
1369:   len               = strlen(name);
1370:   options->names[n] = (char *)malloc((len + 1) * sizeof(char));
1371:   PetscCheck(options->names[n], PETSC_COMM_SELF, PETSC_ERR_MEM, "Failed to allocate option name");
1372:   strcpy(options->names[n], name);

1374: setvalue:
1375:   /* set new value */
1376:   if (options->values[n]) free(options->values[n]);
1377:   len = value ? strlen(value) : 0;
1378:   if (len) {
1379:     options->values[n] = (char *)malloc((len + 1) * sizeof(char));
1380:     if (!options->values[n]) return PETSC_ERR_MEM;
1381:     strcpy(options->values[n], value);
1382:     options->values[n][len] = '\0';
1383:   } else {
1384:     options->values[n] = NULL;
1385:   }
1386:   options->source[n] = source;

1388:   /* handle -help so that it can be set from anywhere */
1389:   if (!PetscOptNameCmp(name, "help")) {
1390:     options->help       = PETSC_TRUE;
1391:     options->help_intro = (value && !PetscOptNameCmp(value, "intro")) ? PETSC_TRUE : PETSC_FALSE;
1392:     options->used[n]    = PETSC_TRUE;
1393:   }

1395:   PetscCall(PetscOptionsMonitor(options, name, value ? value : "", source));
1396:   if (pos) *pos = n;
1397:   PetscFunctionReturn(PETSC_SUCCESS);
1398: }

1400: /*@
1401:   PetscOptionsClearValue - Clears an option name-value pair in the options
1402:   database, overriding whatever is already present.

1404:   Logically Collective

1406:   Input Parameters:
1407: + options - options database, use `NULL` for the default global database
1408: - name    - name of option, this SHOULD have the - prepended

1410:   Level: intermediate

1412:   Note:
1413:   The collectivity of this routine is complex; only the MPI processes that call this routine will
1414:   have the affect of these options. If some processes that create objects call this routine and others do
1415:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1416:   on different ranks.

1418:   Developer Note:
1419:   Uses `free()` directly because the options have been set with `malloc()`

1421: .seealso: `PetscOptionsInsert()`
1422: @*/
1423: PetscErrorCode PetscOptionsClearValue(PetscOptions options, const char name[])
1424: {
1425:   int    N, n, i;
1426:   char **names;

1428:   PetscFunctionBegin;
1429:   options = options ? options : defaultoptions;
1430:   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);
1431:   if (!PetscOptNameCmp(name, "-help")) options->help = options->help_intro = PETSC_FALSE;

1433:   name++; /* skip starting dash */

1435:   /* slow search */
1436:   N = n = options->N;
1437:   names = options->names;
1438:   for (i = 0; i < N; i++) {
1439:     int result = PetscOptNameCmp(names[i], name);
1440:     if (!result) {
1441:       n = i;
1442:       break;
1443:     } else if (result > 0) {
1444:       n = N;
1445:       break;
1446:     }
1447:   }
1448:   if (n == N) PetscFunctionReturn(PETSC_SUCCESS); /* it was not present */

1450:   /* remove name and value */
1451:   if (options->names[n]) free(options->names[n]);
1452:   if (options->values[n]) free(options->values[n]);
1453:   /* shift remaining values down 1 */
1454:   for (i = n; i < N - 1; i++) {
1455:     options->names[i]  = options->names[i + 1];
1456:     options->values[i] = options->values[i + 1];
1457:     options->used[i]   = options->used[i + 1];
1458:     options->source[i] = options->source[i + 1];
1459:   }
1460:   options->N--;

1462:   /* destroy hash table */
1463:   kh_destroy(HO, options->ht);
1464:   options->ht = NULL;

1466:   PetscCall(PetscOptionsMonitor(options, name, NULL, PETSC_OPT_CODE));
1467:   PetscFunctionReturn(PETSC_SUCCESS);
1468: }

1470: /*@C
1471:   PetscOptionsFindPair - Gets an option name-value pair from the options database.

1473:   Not Collective

1475:   Input Parameters:
1476: + options - options database, use `NULL` for the default global database
1477: . pre     - the string to prepend to the name or `NULL`, this SHOULD NOT have the "-" prepended
1478: - name    - name of option, this SHOULD have the "-" prepended

1480:   Output Parameters:
1481: + value - the option value (optional, not used for all options)
1482: - set   - whether the option is set (optional)

1484:   Level: developer

1486:   Note:
1487:   Each process may find different values or no value depending on how options were inserted into the database

1489: .seealso: `PetscOptionsSetValue()`, `PetscOptionsClearValue()`
1490: @*/
1491: PetscErrorCode PetscOptionsFindPair(PetscOptions options, const char pre[], const char name[], const char *value[], PetscBool *set)
1492: {
1493:   char      buf[PETSC_MAX_OPTION_NAME];
1494:   PetscBool usehashtable = PETSC_TRUE;
1495:   PetscBool matchnumbers = PETSC_TRUE;

1497:   PetscFunctionBegin;
1498:   options = options ? options : defaultoptions;
1499:   PetscCheck(!pre || !PetscUnlikely(pre[0] == '-'), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1500:   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);

1502:   name++; /* skip starting dash */

1504:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1505:   if (pre && pre[0]) {
1506:     char *ptr = buf;
1507:     if (name[0] == '-') {
1508:       *ptr++ = '-';
1509:       name++;
1510:     }
1511:     PetscCall(PetscStrncpy(ptr, pre, buf + sizeof(buf) - ptr));
1512:     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
1513:     name = buf;
1514:   }

1516:   if (PetscDefined(USE_DEBUG)) {
1517:     PetscBool valid;
1518:     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
1519:     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
1520:     PetscCall(PetscOptionsValidKey(key, &valid));
1521:     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1522:   }

1524:   if (!options->ht && usehashtable) {
1525:     int          i, ret;
1526:     khiter_t     it;
1527:     khash_t(HO) *ht;
1528:     ht = kh_init(HO);
1529:     PetscCheck(ht, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1530:     ret = kh_resize(HO, ht, options->N * 2); /* twice the required size to reduce risk of collisions */
1531:     PetscCheck(!ret, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1532:     for (i = 0; i < options->N; i++) {
1533:       it = kh_put(HO, ht, options->names[i], &ret);
1534:       PetscCheck(ret == 1, PETSC_COMM_SELF, PETSC_ERR_MEM, "Hash table allocation failed");
1535:       kh_val(ht, it) = i;
1536:     }
1537:     options->ht = ht;
1538:   }

1540:   if (usehashtable) { /* fast search */
1541:     khash_t(HO) *ht = options->ht;
1542:     khiter_t     it = kh_get(HO, ht, name);
1543:     if (it != kh_end(ht)) {
1544:       int i            = kh_val(ht, it);
1545:       options->used[i] = PETSC_TRUE;
1546:       if (value) *value = options->values[i];
1547:       if (set) *set = PETSC_TRUE;
1548:       PetscFunctionReturn(PETSC_SUCCESS);
1549:     }
1550:   } else { /* slow search */
1551:     int i, N = options->N;
1552:     for (i = 0; i < N; i++) {
1553:       int result = PetscOptNameCmp(options->names[i], name);
1554:       if (!result) {
1555:         options->used[i] = PETSC_TRUE;
1556:         if (value) *value = options->values[i];
1557:         if (set) *set = PETSC_TRUE;
1558:         PetscFunctionReturn(PETSC_SUCCESS);
1559:       } else if (result > 0) {
1560:         break;
1561:       }
1562:     }
1563:   }

1565:   /*
1566:    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
1567:    Maybe this special lookup mode should be enabled on request with a push/pop API.
1568:    The feature of matching _%d_ used sparingly in the codebase.
1569:    */
1570:   if (matchnumbers) {
1571:     int i, j, cnt = 0, locs[16], loce[16];
1572:     /* determine the location and number of all _%d_ in the key */
1573:     for (i = 0; name[i]; i++) {
1574:       if (name[i] == '_') {
1575:         for (j = i + 1; name[j]; j++) {
1576:           if (name[j] >= '0' && name[j] <= '9') continue;
1577:           if (name[j] == '_' && j > i + 1) { /* found a number */
1578:             locs[cnt]   = i + 1;
1579:             loce[cnt++] = j + 1;
1580:           }
1581:           i = j - 1;
1582:           break;
1583:         }
1584:       }
1585:     }
1586:     for (i = 0; i < cnt; i++) {
1587:       PetscBool found;
1588:       char      opt[PETSC_MAX_OPTION_NAME + 1] = "-", tmp[PETSC_MAX_OPTION_NAME];
1589:       PetscCall(PetscStrncpy(tmp, name, PetscMin((size_t)(locs[i] + 1), sizeof(tmp))));
1590:       PetscCall(PetscStrlcat(opt, tmp, sizeof(opt)));
1591:       PetscCall(PetscStrlcat(opt, name + loce[i], sizeof(opt)));
1592:       PetscCall(PetscOptionsFindPair(options, NULL, opt, value, &found));
1593:       if (found) {
1594:         if (set) *set = PETSC_TRUE;
1595:         PetscFunctionReturn(PETSC_SUCCESS);
1596:       }
1597:     }
1598:   }

1600:   if (set) *set = PETSC_FALSE;
1601:   PetscFunctionReturn(PETSC_SUCCESS);
1602: }

1604: /* Check whether any option begins with pre+name */
1605: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options, const char pre[], const char name[], const char *option[], const char *value[], PetscBool *set)
1606: {
1607:   char buf[PETSC_MAX_OPTION_NAME];
1608:   int  numCnt = 0, locs[16], loce[16];

1610:   PetscFunctionBegin;
1611:   options = options ? options : defaultoptions;
1612:   PetscCheck(!pre || pre[0] != '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Prefix cannot begin with '-': Instead %s", pre);
1613:   PetscCheck(name[0] == '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Name must begin with '-': Instead %s", name);

1615:   name++; /* skip starting dash */

1617:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1618:   if (pre && pre[0]) {
1619:     char *ptr = buf;
1620:     if (name[0] == '-') {
1621:       *ptr++ = '-';
1622:       name++;
1623:     }
1624:     PetscCall(PetscStrncpy(ptr, pre, sizeof(buf) - ((ptr == buf) ? 0 : 1)));
1625:     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
1626:     name = buf;
1627:   }

1629:   if (PetscDefined(USE_DEBUG)) {
1630:     PetscBool valid;
1631:     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
1632:     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
1633:     PetscCall(PetscOptionsValidKey(key, &valid));
1634:     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1635:   }

1637:   /* determine the location and number of all _%d_ in the key */
1638:   {
1639:     int i, j;
1640:     for (i = 0; name[i]; i++) {
1641:       if (name[i] == '_') {
1642:         for (j = i + 1; name[j]; j++) {
1643:           if (name[j] >= '0' && name[j] <= '9') continue;
1644:           if (name[j] == '_' && j > i + 1) { /* found a number */
1645:             locs[numCnt]   = i + 1;
1646:             loce[numCnt++] = j + 1;
1647:           }
1648:           i = j - 1;
1649:           break;
1650:         }
1651:       }
1652:     }
1653:   }

1655:   /* slow search */
1656:   for (int c = -1; c < numCnt; ++c) {
1657:     char   opt[PETSC_MAX_OPTION_NAME + 2] = "";
1658:     size_t len;

1660:     if (c < 0) {
1661:       PetscCall(PetscStrncpy(opt, name, sizeof(opt)));
1662:     } else {
1663:       PetscCall(PetscStrncpy(opt, name, PetscMin((size_t)(locs[c] + 1), sizeof(opt))));
1664:       PetscCall(PetscStrlcat(opt, name + loce[c], sizeof(opt) - 1));
1665:     }
1666:     PetscCall(PetscStrlen(opt, &len));
1667:     for (int i = 0; i < options->N; i++) {
1668:       PetscBool match;

1670:       PetscCall(PetscStrncmp(options->names[i], opt, len, &match));
1671:       if (match) {
1672:         options->used[i] = PETSC_TRUE;
1673:         if (option) *option = options->names[i];
1674:         if (value) *value = options->values[i];
1675:         if (set) *set = PETSC_TRUE;
1676:         PetscFunctionReturn(PETSC_SUCCESS);
1677:       }
1678:     }
1679:   }

1681:   if (set) *set = PETSC_FALSE;
1682:   PetscFunctionReturn(PETSC_SUCCESS);
1683: }

1685: /*@
1686:   PetscOptionsReject - Generates an error if a certain option is given.

1688:   Not Collective

1690:   Input Parameters:
1691: + options - options database, use `NULL` for default global database
1692: . pre     - the option prefix (may be `NULL`)
1693: . name    - the option name one is seeking
1694: - mess    - error message (may be `NULL`)

1696:   Level: advanced

1698: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`,
1699:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1700:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1701:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1702:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1703:           `PetscOptionsFList()`, `PetscOptionsEList()`
1704: @*/
1705: PetscErrorCode PetscOptionsReject(PetscOptions options, const char pre[], const char name[], const char mess[])
1706: {
1707:   PetscBool flag = PETSC_FALSE;

1709:   PetscFunctionBegin;
1710:   PetscCall(PetscOptionsHasName(options, pre, name, &flag));
1711:   if (flag) {
1712:     PetscCheck(!mess || !mess[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s with %s", pre ? pre : "", name + 1, mess);
1713:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s", pre ? pre : "", name + 1);
1714:   }
1715:   PetscFunctionReturn(PETSC_SUCCESS);
1716: }

1718: /*@
1719:   PetscOptionsHasHelp - Determines whether the "-help" option is in the database.

1721:   Not Collective

1723:   Input Parameter:
1724: . options - options database, use `NULL` for default global database

1726:   Output Parameter:
1727: . set - `PETSC_TRUE` if found else `PETSC_FALSE`.

1729:   Level: advanced

1731: .seealso: `PetscOptionsHasName()`
1732: @*/
1733: PetscErrorCode PetscOptionsHasHelp(PetscOptions options, PetscBool *set)
1734: {
1735:   PetscFunctionBegin;
1736:   PetscAssertPointer(set, 2);
1737:   options = options ? options : defaultoptions;
1738:   *set    = options->help;
1739:   PetscFunctionReturn(PETSC_SUCCESS);
1740: }

1742: PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options, PetscBool *set)
1743: {
1744:   PetscFunctionBegin;
1745:   PetscAssertPointer(set, 2);
1746:   options = options ? options : defaultoptions;
1747:   *set    = options->help_intro;
1748:   PetscFunctionReturn(PETSC_SUCCESS);
1749: }

1751: /*@
1752:   PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or Boolean, even
1753:   if its value is set to false.

1755:   Not Collective

1757:   Input Parameters:
1758: + options - options database, use `NULL` for default global database
1759: . pre     - string to prepend to the name or `NULL`
1760: - name    - the option one is seeking

1762:   Output Parameter:
1763: . set - `PETSC_TRUE` if found else `PETSC_FALSE`.

1765:   Level: beginner

1767:   Note:
1768:   In many cases you probably want to use `PetscOptionsGetBool()` instead of calling this, to allowing toggling values.

1770: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
1771:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1772:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1773:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1774:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1775:           `PetscOptionsFList()`, `PetscOptionsEList()`
1776: @*/
1777: PetscErrorCode PetscOptionsHasName(PetscOptions options, const char pre[], const char name[], PetscBool *set)
1778: {
1779:   const char *value;
1780:   PetscBool   flag;

1782:   PetscFunctionBegin;
1783:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
1784:   if (set) *set = flag;
1785:   PetscFunctionReturn(PETSC_SUCCESS);
1786: }

1788: /*@C
1789:   PetscOptionsGetAll - Lists all the options the program was run with in a single string.

1791:   Not Collective

1793:   Input Parameter:
1794: . options - the options database, use `NULL` for the default global database

1796:   Output Parameter:
1797: . copts - pointer where string pointer is stored

1799:   Level: advanced

1801:   Notes:
1802:   The array and each entry in the array should be freed with `PetscFree()`

1804:   Each process may have different values depending on how the options were inserted into the database

1806: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsView()`, `PetscOptionsPush()`, `PetscOptionsPop()`
1807: @*/
1808: PetscErrorCode PetscOptionsGetAll(PetscOptions options, char *copts[])
1809: {
1810:   PetscInt i;
1811:   size_t   len = 1, lent = 0;
1812:   char    *coptions = NULL;

1814:   PetscFunctionBegin;
1815:   PetscAssertPointer(copts, 2);
1816:   options = options ? options : defaultoptions;
1817:   /* count the length of the required string */
1818:   for (i = 0; i < options->N; i++) {
1819:     PetscCall(PetscStrlen(options->names[i], &lent));
1820:     len += 2 + lent;
1821:     if (options->values[i]) {
1822:       PetscCall(PetscStrlen(options->values[i], &lent));
1823:       len += 1 + lent;
1824:     }
1825:   }
1826:   PetscCall(PetscMalloc1(len, &coptions));
1827:   coptions[0] = 0;
1828:   for (i = 0; i < options->N; i++) {
1829:     PetscCall(PetscStrlcat(coptions, "-", len));
1830:     PetscCall(PetscStrlcat(coptions, options->names[i], len));
1831:     PetscCall(PetscStrlcat(coptions, " ", len));
1832:     if (options->values[i]) {
1833:       PetscCall(PetscStrlcat(coptions, options->values[i], len));
1834:       PetscCall(PetscStrlcat(coptions, " ", len));
1835:     }
1836:   }
1837:   *copts = coptions;
1838:   PetscFunctionReturn(PETSC_SUCCESS);
1839: }

1841: /*@
1842:   PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database

1844:   Not Collective

1846:   Input Parameters:
1847: + options - options database, use `NULL` for default global database
1848: - name    - string name of option

1850:   Output Parameter:
1851: . used - `PETSC_TRUE` if the option was used, otherwise false, including if option was not found in options database

1853:   Level: advanced

1855:   Note:
1856:   The value returned may be different on each process and depends on which options have been processed
1857:   on the given process

1859: .seealso: `PetscOptionsView()`, `PetscOptionsLeft()`, `PetscOptionsAllUsed()`
1860: @*/
1861: PetscErrorCode PetscOptionsUsed(PetscOptions options, const char *name, PetscBool *used)
1862: {
1863:   PetscInt i;

1865:   PetscFunctionBegin;
1866:   PetscAssertPointer(name, 2);
1867:   PetscAssertPointer(used, 3);
1868:   options = options ? options : defaultoptions;
1869:   *used   = PETSC_FALSE;
1870:   for (i = 0; i < options->N; i++) {
1871:     PetscCall(PetscStrcasecmp(options->names[i], name, used));
1872:     if (*used) {
1873:       *used = options->used[i];
1874:       break;
1875:     }
1876:   }
1877:   PetscFunctionReturn(PETSC_SUCCESS);
1878: }

1880: /*@
1881:   PetscOptionsAllUsed - Returns a count of the number of options in the
1882:   database that have never been selected.

1884:   Not Collective

1886:   Input Parameter:
1887: . options - options database, use `NULL` for default global database

1889:   Output Parameter:
1890: . N - count of options not used

1892:   Level: advanced

1894:   Note:
1895:   The value returned may be different on each process and depends on which options have been processed
1896:   on the given process

1898: .seealso: `PetscOptionsView()`
1899: @*/
1900: PetscErrorCode PetscOptionsAllUsed(PetscOptions options, PetscInt *N)
1901: {
1902:   PetscInt i, n = 0;

1904:   PetscFunctionBegin;
1905:   PetscAssertPointer(N, 2);
1906:   options = options ? options : defaultoptions;
1907:   for (i = 0; i < options->N; i++) {
1908:     if (!options->used[i]) n++;
1909:   }
1910:   *N = n;
1911:   PetscFunctionReturn(PETSC_SUCCESS);
1912: }

1914: /*@
1915:   PetscOptionsLeft - Prints to screen any options that were set and never used.

1917:   Not Collective

1919:   Input Parameter:
1920: . options - options database; use `NULL` for default global database

1922:   Options Database Key:
1923: . -options_left - activates `PetscOptionsAllUsed()` within `PetscFinalize()`

1925:   Level: advanced

1927:   Notes:
1928:   This is rarely used directly, it is called by `PetscFinalize()` in debug more or if -options_left
1929:   is passed otherwise to help users determine possible mistakes in their usage of options. This
1930:   only prints values on process zero of `PETSC_COMM_WORLD`.

1932:   Other processes depending the objects
1933:   used may have different options that are left unused.

1935: .seealso: `PetscOptionsAllUsed()`
1936: @*/
1937: PetscErrorCode PetscOptionsLeft(PetscOptions options)
1938: {
1939:   PetscInt     i;
1940:   PetscInt     cnt = 0;
1941:   PetscOptions toptions;

1943:   PetscFunctionBegin;
1944:   toptions = options ? options : defaultoptions;
1945:   for (i = 0; i < toptions->N; i++) {
1946:     if (!toptions->used[i]) {
1947:       if (PetscCIOption(toptions->names[i])) continue;
1948:       if (toptions->values[i]) {
1949:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s value: %s source: %s\n", toptions->names[i], toptions->values[i], PetscOptionSources[toptions->source[i]]));
1950:       } else {
1951:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s (no value) source: %s\n", toptions->names[i], PetscOptionSources[toptions->source[i]]));
1952:       }
1953:     }
1954:   }
1955:   if (!options) {
1956:     toptions = defaultoptions;
1957:     while (toptions->previous) {
1958:       cnt++;
1959:       toptions = toptions->previous;
1960:     }
1961:     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));
1962:   }
1963:   PetscFunctionReturn(PETSC_SUCCESS);
1964: }

1966: /*@C
1967:   PetscOptionsLeftGet - Returns all options that were set and never used.

1969:   Not Collective

1971:   Input Parameter:
1972: . options - options database, use `NULL` for default global database

1974:   Output Parameters:
1975: + N      - count of options not used
1976: . names  - names of options not used
1977: - values - values of options not used

1979:   Level: advanced

1981:   Notes:
1982:   Users should call `PetscOptionsLeftRestore()` to free the memory allocated in this routine

1984:   The value returned may be different on each process and depends on which options have been processed
1985:   on the given process

1987: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`
1988: @*/
1989: PetscErrorCode PetscOptionsLeftGet(PetscOptions options, PetscInt *N, char **names[], char **values[])
1990: {
1991:   PetscInt i, n;

1993:   PetscFunctionBegin;
1994:   if (N) PetscAssertPointer(N, 2);
1995:   if (names) PetscAssertPointer(names, 3);
1996:   if (values) PetscAssertPointer(values, 4);
1997:   options = options ? options : defaultoptions;

1999:   /* The number of unused PETSc options */
2000:   n = 0;
2001:   for (i = 0; i < options->N; i++) {
2002:     if (PetscCIOption(options->names[i])) continue;
2003:     if (!options->used[i]) n++;
2004:   }
2005:   if (N) *N = n;
2006:   if (names) PetscCall(PetscMalloc1(n, names));
2007:   if (values) PetscCall(PetscMalloc1(n, values));

2009:   n = 0;
2010:   if (names || values) {
2011:     for (i = 0; i < options->N; i++) {
2012:       if (!options->used[i]) {
2013:         if (PetscCIOption(options->names[i])) continue;
2014:         if (names) (*names)[n] = options->names[i];
2015:         if (values) (*values)[n] = options->values[i];
2016:         n++;
2017:       }
2018:     }
2019:   }
2020:   PetscFunctionReturn(PETSC_SUCCESS);
2021: }

2023: /*@C
2024:   PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using `PetscOptionsLeftGet()`.

2026:   Not Collective

2028:   Input Parameters:
2029: + options - options database, use `NULL` for default global database
2030: . N       - count of options not used
2031: . names   - names of options not used
2032: - values  - values of options not used

2034:   Level: advanced

2036:   Notes:
2037:   The user should pass the same pointer to `N` as they did when calling `PetscOptionsLeftGet()`

2039: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`, `PetscOptionsLeftGet()`
2040: @*/
2041: PetscErrorCode PetscOptionsLeftRestore(PetscOptions options, PetscInt *N, char **names[], char **values[])
2042: {
2043:   PetscFunctionBegin;
2044:   (void)options;
2045:   if (N) PetscAssertPointer(N, 2);
2046:   if (names) PetscAssertPointer(names, 3);
2047:   if (values) PetscAssertPointer(values, 4);
2048:   if (N) *N = 0;
2049:   if (names) PetscCall(PetscFree(*names));
2050:   if (values) PetscCall(PetscFree(*values));
2051:   PetscFunctionReturn(PETSC_SUCCESS);
2052: }

2054: /*@C
2055:   PetscOptionsMonitorDefault - Print all options set value events using the supplied `PetscViewer`.

2057:   Logically Collective

2059:   Input Parameters:
2060: + name   - option name string
2061: . value  - option value string
2062: . source - The source for the option
2063: - ctx    - a `PETSCVIEWERASCII` or `NULL`

2065:   Level: intermediate

2067:   Notes:
2068:   If ctx is `NULL`, `PetscPrintf()` is used.
2069:   The first MPI process in the `PetscViewer` viewer actually prints the values, other
2070:   processes may have different values set

2072:   If `PetscCIEnabled` then do not print the test harness options

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);

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 (source: %s)\n", name, value, PetscOptionSources[source]));
2098:     }
2099:   }
2100:   PetscFunctionReturn(PETSC_SUCCESS);
2101: }

2103: /*@C
2104:   PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
2105:   modified the PETSc options database.

2107:   Not Collective

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 monitor routine (use `NULL` if
2112:                    no context is desired)
2113: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

2115:   Calling sequence of `monitor`:
2116: + name   - option name string
2117: . value  - option value string, a value of `NULL` indicates the option is being removed from the database. A value
2118:            of "" indicates the option is in the database but has no value.
2119: . source - option source
2120: - mctx   - optional monitoring context, as set by `PetscOptionsMonitorSet()`

2122:   Calling sequence of `monitordestroy`:
2123: . mctx - [optional] pointer to context to destroy with

2125:   Options Database Keys:
2126: + -options_monitor <viewer> - turn on default monitoring
2127: - -options_monitor_cancel   - turn off any option monitors except the default monitor obtained with `-options_monitor`

2129:   Level: intermediate

2131:   Notes:
2132:   See `PetscInitialize()` for options related to option database monitoring.

2134:   The default is to do no monitoring.  To print the name and value of options
2135:   being inserted into the database, use `PetscOptionsMonitorDefault()` as the monitoring routine,
2136:   with a `NULL` monitoring context. Or use the option `-options_monitor` <viewer>.

2138:   Several different monitoring routines may be set by calling
2139:   `PetscOptionsMonitorSet()` multiple times; all will be called in the
2140:   order in which they were set.

2142: .seealso: `PetscOptionsMonitorDefault()`, `PetscInitialize()`
2143: @*/
2144: PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], PetscOptionSource source, void *mctx), void *mctx, PetscErrorCode (*monitordestroy)(void **mctx))
2145: {
2146:   PetscOptions options = defaultoptions;

2148:   PetscFunctionBegin;
2149:   if (options->monitorCancel) PetscFunctionReturn(PETSC_SUCCESS);
2150:   PetscCheck(options->numbermonitors < MAXOPTIONSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many PetscOptions monitors set");
2151:   options->monitor[options->numbermonitors]          = monitor;
2152:   options->monitordestroy[options->numbermonitors]   = monitordestroy;
2153:   options->monitorcontext[options->numbermonitors++] = (void *)mctx;
2154:   PetscFunctionReturn(PETSC_SUCCESS);
2155: }

2157: /*
2158:    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
2159:      Empty string is considered as true.
2160: */
2161: PetscErrorCode PetscOptionsStringToBool(const char value[], PetscBool *a)
2162: {
2163:   PetscBool istrue, isfalse;
2164:   size_t    len;

2166:   PetscFunctionBegin;
2167:   /* PetscStrlen() returns 0 for NULL or "" */
2168:   PetscCall(PetscStrlen(value, &len));
2169:   if (!len) {
2170:     *a = PETSC_TRUE;
2171:     PetscFunctionReturn(PETSC_SUCCESS);
2172:   }
2173:   PetscCall(PetscStrcasecmp(value, "TRUE", &istrue));
2174:   if (istrue) {
2175:     *a = PETSC_TRUE;
2176:     PetscFunctionReturn(PETSC_SUCCESS);
2177:   }
2178:   PetscCall(PetscStrcasecmp(value, "YES", &istrue));
2179:   if (istrue) {
2180:     *a = PETSC_TRUE;
2181:     PetscFunctionReturn(PETSC_SUCCESS);
2182:   }
2183:   PetscCall(PetscStrcasecmp(value, "1", &istrue));
2184:   if (istrue) {
2185:     *a = PETSC_TRUE;
2186:     PetscFunctionReturn(PETSC_SUCCESS);
2187:   }
2188:   PetscCall(PetscStrcasecmp(value, "on", &istrue));
2189:   if (istrue) {
2190:     *a = PETSC_TRUE;
2191:     PetscFunctionReturn(PETSC_SUCCESS);
2192:   }
2193:   PetscCall(PetscStrcasecmp(value, "FALSE", &isfalse));
2194:   if (isfalse) {
2195:     *a = PETSC_FALSE;
2196:     PetscFunctionReturn(PETSC_SUCCESS);
2197:   }
2198:   PetscCall(PetscStrcasecmp(value, "NO", &isfalse));
2199:   if (isfalse) {
2200:     *a = PETSC_FALSE;
2201:     PetscFunctionReturn(PETSC_SUCCESS);
2202:   }
2203:   PetscCall(PetscStrcasecmp(value, "0", &isfalse));
2204:   if (isfalse) {
2205:     *a = PETSC_FALSE;
2206:     PetscFunctionReturn(PETSC_SUCCESS);
2207:   }
2208:   PetscCall(PetscStrcasecmp(value, "off", &isfalse));
2209:   if (isfalse) {
2210:     *a = PETSC_FALSE;
2211:     PetscFunctionReturn(PETSC_SUCCESS);
2212:   }
2213:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown logical value: %s", value);
2214: }

2216: /*
2217:    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2218: */
2219: PetscErrorCode PetscOptionsStringToInt(const char name[], PetscInt *a)
2220: {
2221:   size_t    len;
2222:   PetscBool decide, tdefault, mouse, unlimited;

2224:   PetscFunctionBegin;
2225:   PetscCall(PetscStrlen(name, &len));
2226:   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");

2228:   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &tdefault));
2229:   if (!tdefault) PetscCall(PetscStrcasecmp(name, "DEFAULT", &tdefault));
2230:   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &decide));
2231:   if (!decide) PetscCall(PetscStrcasecmp(name, "DECIDE", &decide));
2232:   if (!decide) PetscCall(PetscStrcasecmp(name, "PETSC_DETERMINE", &decide));
2233:   if (!decide) PetscCall(PetscStrcasecmp(name, "DETERMINE", &decide));
2234:   PetscCall(PetscStrcasecmp(name, "PETSC_UNLIMITED", &unlimited));
2235:   if (!unlimited) PetscCall(PetscStrcasecmp(name, "UNLIMITED", &unlimited));
2236:   PetscCall(PetscStrcasecmp(name, "mouse", &mouse));

2238:   if (tdefault) *a = PETSC_DEFAULT;
2239:   else if (decide) *a = PETSC_DECIDE;
2240:   else if (unlimited) *a = PETSC_UNLIMITED;
2241:   else if (mouse) *a = -1;
2242:   else {
2243:     char *endptr;
2244:     long  strtolval;

2246:     strtolval = strtol(name, &endptr, 10);
2247:     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);

2249: #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2250:     (void)strtolval;
2251:     *a = atoll(name);
2252: #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2253:     (void)strtolval;
2254:     *a = _atoi64(name);
2255: #else
2256:     *a = (PetscInt)strtolval;
2257: #endif
2258:   }
2259:   PetscFunctionReturn(PETSC_SUCCESS);
2260: }

2262: #if defined(PETSC_USE_REAL___FLOAT128)
2263:   #include <quadmath.h>
2264: #endif

2266: static PetscErrorCode PetscStrtod(const char name[], PetscReal *a, char **endptr)
2267: {
2268:   PetscFunctionBegin;
2269: #if defined(PETSC_USE_REAL___FLOAT128)
2270:   *a = strtoflt128(name, endptr);
2271: #else
2272:   *a = (PetscReal)strtod(name, endptr);
2273: #endif
2274:   PetscFunctionReturn(PETSC_SUCCESS);
2275: }

2277: static PetscErrorCode PetscStrtoz(const char name[], PetscScalar *a, char **endptr, PetscBool *isImaginary)
2278: {
2279:   PetscBool hasi = PETSC_FALSE;
2280:   char     *ptr;
2281:   PetscReal strtoval;

2283:   PetscFunctionBegin;
2284:   PetscCall(PetscStrtod(name, &strtoval, &ptr));
2285:   if (ptr == name) {
2286:     strtoval = 1.;
2287:     hasi     = PETSC_TRUE;
2288:     if (name[0] == 'i') {
2289:       ptr++;
2290:     } else if (name[0] == '+' && name[1] == 'i') {
2291:       ptr += 2;
2292:     } else if (name[0] == '-' && name[1] == 'i') {
2293:       strtoval = -1.;
2294:       ptr += 2;
2295:     }
2296:   } else if (*ptr == 'i') {
2297:     hasi = PETSC_TRUE;
2298:     ptr++;
2299:   }
2300:   *endptr      = ptr;
2301:   *isImaginary = hasi;
2302:   if (hasi) {
2303: #if !defined(PETSC_USE_COMPLEX)
2304:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s contains imaginary but complex not supported ", name);
2305: #else
2306:     *a = PetscCMPLX(0., strtoval);
2307: #endif
2308:   } else {
2309:     *a = strtoval;
2310:   }
2311:   PetscFunctionReturn(PETSC_SUCCESS);
2312: }

2314: /*
2315:    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2316: */
2317: PetscErrorCode PetscOptionsStringToReal(const char name[], PetscReal *a)
2318: {
2319:   size_t    len;
2320:   PetscBool match;
2321:   char     *endptr;

2323:   PetscFunctionBegin;
2324:   PetscCall(PetscStrlen(name, &len));
2325:   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "String of length zero has no numerical value");

2327:   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &match));
2328:   if (!match) PetscCall(PetscStrcasecmp(name, "DEFAULT", &match));
2329:   if (match) {
2330:     *a = PETSC_DEFAULT;
2331:     PetscFunctionReturn(PETSC_SUCCESS);
2332:   }

2334:   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &match));
2335:   if (!match) PetscCall(PetscStrcasecmp(name, "DECIDE", &match));
2336:   if (match) {
2337:     *a = PETSC_DECIDE;
2338:     PetscFunctionReturn(PETSC_SUCCESS);
2339:   }

2341:   PetscCall(PetscStrcasecmp(name, "PETSC_DETERMINE", &match));
2342:   if (!match) PetscCall(PetscStrcasecmp(name, "DETERMINE", &match));
2343:   if (match) {
2344:     *a = PETSC_DETERMINE;
2345:     PetscFunctionReturn(PETSC_SUCCESS);
2346:   }

2348:   PetscCall(PetscStrcasecmp(name, "PETSC_UNLIMITED", &match));
2349:   if (!match) PetscCall(PetscStrcasecmp(name, "UNLIMITED", &match));
2350:   if (match) {
2351:     *a = PETSC_UNLIMITED;
2352:     PetscFunctionReturn(PETSC_SUCCESS);
2353:   }

2355:   PetscCall(PetscStrtod(name, a, &endptr));
2356:   PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value", name);
2357:   PetscFunctionReturn(PETSC_SUCCESS);
2358: }

2360: PetscErrorCode PetscOptionsStringToScalar(const char name[], PetscScalar *a)
2361: {
2362:   PetscBool   imag1;
2363:   size_t      len;
2364:   PetscScalar val = 0.;
2365:   char       *ptr = NULL;

2367:   PetscFunctionBegin;
2368:   PetscCall(PetscStrlen(name, &len));
2369:   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
2370:   PetscCall(PetscStrtoz(name, &val, &ptr, &imag1));
2371: #if defined(PETSC_USE_COMPLEX)
2372:   if ((size_t)(ptr - name) < len) {
2373:     PetscBool   imag2;
2374:     PetscScalar val2;

2376:     PetscCall(PetscStrtoz(ptr, &val2, &ptr, &imag2));
2377:     if (imag1) PetscCheck(imag2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s: must specify imaginary component second", name);
2378:     val = PetscCMPLX(PetscRealPart(val), PetscImaginaryPart(val2));
2379:   }
2380: #endif
2381:   PetscCheck((size_t)(ptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value ", name);
2382:   *a = val;
2383:   PetscFunctionReturn(PETSC_SUCCESS);
2384: }

2386: /*@C
2387:   PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2388:   option in the database.

2390:   Not Collective

2392:   Input Parameters:
2393: + options - options database, use `NULL` for default global database
2394: . pre     - the string to prepend to the name or `NULL`
2395: - name    - the option one is seeking

2397:   Output Parameters:
2398: + ivalue - the logical value to return
2399: - set    - `PETSC_TRUE`  if found, else `PETSC_FALSE`

2401:   Level: beginner

2403:   Notes:
2404:   TRUE, true, YES, yes, ON, on, nostring, and 1 all translate to `PETSC_TRUE`
2405:   FALSE, false, NO, no, OFF, off and 0 all translate to `PETSC_FALSE`

2407:   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`
2408:   is equivalent to `-requested_bool true`

2410:   If the user does not supply the option at all `ivalue` is NOT changed. Thus
2411:   you should ALWAYS initialize `ivalue` if you access it without first checking that the `set` flag is true.

2413: .seealso: `PetscOptionsGetBool3()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2414:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2415:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2416:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2417:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2418:           `PetscOptionsFList()`, `PetscOptionsEList()`
2419: @*/
2420: PetscErrorCode PetscOptionsGetBool(PetscOptions options, const char pre[], const char name[], PetscBool *ivalue, PetscBool *set)
2421: {
2422:   const char *value;
2423:   PetscBool   flag;

2425:   PetscFunctionBegin;
2426:   PetscAssertPointer(name, 3);
2427:   if (ivalue) PetscAssertPointer(ivalue, 4);
2428:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2429:   if (flag) {
2430:     if (set) *set = PETSC_TRUE;
2431:     PetscCall(PetscOptionsStringToBool(value, &flag));
2432:     if (ivalue) *ivalue = flag;
2433:   } else {
2434:     if (set) *set = PETSC_FALSE;
2435:   }
2436:   PetscFunctionReturn(PETSC_SUCCESS);
2437: }

2439: /*@C
2440:   PetscOptionsGetBool3 - Gets the ternary logical (true, false or unknown) value for a particular
2441:   option in the database.

2443:   Not Collective

2445:   Input Parameters:
2446: + options - options database, use `NULL` for default global database
2447: . pre     - the string to prepend to the name or `NULL`
2448: - name    - the option one is seeking

2450:   Output Parameters:
2451: + ivalue - the ternary logical value to return
2452: - set    - `PETSC_TRUE`  if found, else `PETSC_FALSE`

2454:   Level: beginner

2456:   Notes:
2457:   TRUE, true, YES, yes, ON, on, nostring and 1 all translate to `PETSC_BOOL3_TRUE`
2458:   FALSE, false, NO, no, OFF, off and 0 all translate to `PETSC_BOOL3_FALSE`
2459:   UNKNOWN, unknown, AUTO and auto all translate to `PETSC_BOOL3_UNKNOWN`

2461:   If the option is given, but no value is provided, then `ivalue` will be set to `PETSC_BOOL3_TRUE` and `set` will be set to `PETSC_TRUE`. That is `-requested_bool3`
2462:   is equivalent to `-requested_bool3 true`

2464:   If the user does not supply the option at all `ivalue` is NOT changed. Thus
2465:   you should ALWAYS initialize `ivalue` if you access it without first checking that the `set` flag is true.

2467: .seealso: `PetscOptionsGetBool()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2468:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2469:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2470:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2471:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2472:           `PetscOptionsFList()`, `PetscOptionsEList()`
2473: @*/
2474: PetscErrorCode PetscOptionsGetBool3(PetscOptions options, const char pre[], const char name[], PetscBool3 *ivalue, PetscBool *set)
2475: {
2476:   const char *value;
2477:   PetscBool   flag;

2479:   PetscFunctionBegin;
2480:   PetscAssertPointer(name, 3);
2481:   if (ivalue) PetscAssertPointer(ivalue, 4);
2482:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2483:   if (flag) { // found the option
2484:     PetscBool isAUTO = PETSC_FALSE, isUNKNOWN = PETSC_FALSE;

2486:     if (set) *set = PETSC_TRUE;
2487:     PetscCall(PetscStrcasecmp("AUTO", value, &isAUTO));                    // auto or AUTO
2488:     if (!isAUTO) PetscCall(PetscStrcasecmp("UNKNOWN", value, &isUNKNOWN)); // unknown or UNKNOWN
2489:     if (isAUTO || isUNKNOWN) {
2490:       if (ivalue) *ivalue = PETSC_BOOL3_UNKNOWN;
2491:     } else { // handle boolean values (if no value is given, it returns true)
2492:       PetscCall(PetscOptionsStringToBool(value, &flag));
2493:       if (ivalue) *ivalue = PetscBoolToBool3(flag);
2494:     }
2495:   } else {
2496:     if (set) *set = PETSC_FALSE;
2497:   }
2498:   PetscFunctionReturn(PETSC_SUCCESS);
2499: }

2501: /*@C
2502:   PetscOptionsGetEList - Puts a list of option values that a single one may be selected from

2504:   Not Collective

2506:   Input Parameters:
2507: + options - options database, use `NULL` for default global database
2508: . pre     - the string to prepend to the name or `NULL`
2509: . opt     - option name
2510: . list    - the possible choices (one of these must be selected, anything else is invalid)
2511: - ntext   - number of choices

2513:   Output Parameters:
2514: + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2515: - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`

2517:   Level: intermediate

2519:   Notes:
2520:   If the user does not supply the option `value` is NOT changed. Thus
2521:   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is true.

2523:   See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList`

2525: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2526:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2527:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2528:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2529:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2530:           `PetscOptionsFList()`, `PetscOptionsEList()`
2531: @*/
2532: PetscErrorCode PetscOptionsGetEList(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscInt ntext, PetscInt *value, PetscBool *set)
2533: {
2534:   size_t    alen, len = 0, tlen = 0;
2535:   char     *svalue;
2536:   PetscBool aset, flg = PETSC_FALSE;
2537:   PetscInt  i;

2539:   PetscFunctionBegin;
2540:   PetscAssertPointer(opt, 3);
2541:   for (i = 0; i < ntext; i++) {
2542:     PetscCall(PetscStrlen(list[i], &alen));
2543:     if (alen > len) len = alen;
2544:     tlen += len + 1;
2545:   }
2546:   len += 5; /* a little extra space for user mistypes */
2547:   PetscCall(PetscMalloc1(len, &svalue));
2548:   PetscCall(PetscOptionsGetString(options, pre, opt, svalue, len, &aset));
2549:   if (aset) {
2550:     PetscCall(PetscEListFind(ntext, list, svalue, value, &flg));
2551:     if (!flg) {
2552:       char *avail;

2554:       PetscCall(PetscMalloc1(tlen, &avail));
2555:       avail[0] = '\0';
2556:       for (i = 0; i < ntext; i++) {
2557:         PetscCall(PetscStrlcat(avail, list[i], tlen));
2558:         PetscCall(PetscStrlcat(avail, " ", tlen));
2559:       }
2560:       PetscCall(PetscStrtolower(avail));
2561:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown option %s for -%s%s. Available options: %s", svalue, pre ? pre : "", opt + 1, avail);
2562:     }
2563:     if (set) *set = PETSC_TRUE;
2564:   } else if (set) *set = PETSC_FALSE;
2565:   PetscCall(PetscFree(svalue));
2566:   PetscFunctionReturn(PETSC_SUCCESS);
2567: }

2569: /*@C
2570:   PetscOptionsGetEnum - Gets the enum value for a particular option in the database.

2572:   Not Collective

2574:   Input Parameters:
2575: + options - options database, use `NULL` for default global database
2576: . pre     - option prefix or `NULL`
2577: . opt     - option name
2578: - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2580:   Output Parameters:
2581: + value - the value to return
2582: - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`

2584:   Level: beginner

2586:   Notes:
2587:   If the user does not supply the option `value` is NOT changed. Thus
2588:   you should ALWAYS initialize `value` if you access it without first checking that the `set` flag is true.

2590:   `list` is usually something like `PCASMTypes` or some other predefined list of enum names

2592: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2593:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2594:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2595:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2596:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2597:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2598:           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2599: @*/
2600: PetscErrorCode PetscOptionsGetEnum(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscEnum *value, PetscBool *set)
2601: {
2602:   PetscInt  ntext = 0, tval;
2603:   PetscBool fset;

2605:   PetscFunctionBegin;
2606:   PetscAssertPointer(opt, 3);
2607:   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");
2608:   PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix");
2609:   ntext -= 3;
2610:   PetscCall(PetscOptionsGetEList(options, pre, opt, list, ntext, &tval, &fset));
2611:   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2612:   if (fset) *value = (PetscEnum)tval;
2613:   if (set) *set = fset;
2614:   PetscFunctionReturn(PETSC_SUCCESS);
2615: }

2617: /*@C
2618:   PetscOptionsGetInt - Gets the integer value for a particular option in the database.

2620:   Not Collective

2622:   Input Parameters:
2623: + options - options database, use `NULL` for default global database
2624: . pre     - the string to prepend to the name or `NULL`
2625: - name    - the option one is seeking

2627:   Output Parameters:
2628: + ivalue - the integer value to return
2629: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2631:   Level: beginner

2633:   Notes:
2634:   If the user does not supply the option `ivalue` is NOT changed. Thus
2635:   you should ALWAYS initialize the `ivalue` if you access it without first checking that the `set` flag is true.

2637:   Accepts the special values `determine`, `decide` and `unlimited`.

2639:   Accepts the deprecated value `default`.

2641: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2642:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2643:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2644:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2645:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2646:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2647:           `PetscOptionsFList()`, `PetscOptionsEList()`
2648: @*/
2649: PetscErrorCode PetscOptionsGetInt(PetscOptions options, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
2650: {
2651:   const char *value;
2652:   PetscBool   flag;

2654:   PetscFunctionBegin;
2655:   PetscAssertPointer(name, 3);
2656:   PetscAssertPointer(ivalue, 4);
2657:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2658:   if (flag) {
2659:     if (!value) {
2660:       if (set) *set = PETSC_FALSE;
2661:     } else {
2662:       if (set) *set = PETSC_TRUE;
2663:       PetscCall(PetscOptionsStringToInt(value, ivalue));
2664:     }
2665:   } else {
2666:     if (set) *set = PETSC_FALSE;
2667:   }
2668:   PetscFunctionReturn(PETSC_SUCCESS);
2669: }

2671: /*@C
2672:   PetscOptionsGetMPIInt - Gets the MPI integer value for a particular option in the database.

2674:   Not Collective

2676:   Input Parameters:
2677: + options - options database, use `NULL` for default global database
2678: . pre     - the string to prepend to the name or `NULL`
2679: - name    - the option one is seeking

2681:   Output Parameters:
2682: + ivalue - the MPI integer value to return
2683: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2685:   Level: beginner

2687:   Notes:
2688:   If the user does not supply the option `ivalue` is NOT changed. Thus
2689:   you should ALWAYS initialize the `ivalue` if you access it without first checking that the `set` flag is true.

2691:   Accepts the special values `determine`, `decide` and `unlimited`.

2693:   Accepts the deprecated value `default`.

2695: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2696:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2697:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2698:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2699:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2700:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2701:           `PetscOptionsFList()`, `PetscOptionsEList()`
2702: @*/
2703: PetscErrorCode PetscOptionsGetMPIInt(PetscOptions options, const char pre[], const char name[], PetscMPIInt *ivalue, PetscBool *set)
2704: {
2705:   PetscInt  value;
2706:   PetscBool flag;

2708:   PetscFunctionBegin;
2709:   PetscCall(PetscOptionsGetInt(options, pre, name, &value, &flag));
2710:   if (flag) PetscCall(PetscMPIIntCast(value, ivalue));
2711:   if (set) *set = flag;
2712:   PetscFunctionReturn(PETSC_SUCCESS);
2713: }

2715: /*@C
2716:   PetscOptionsGetReal - Gets the double precision value for a particular
2717:   option in the database.

2719:   Not Collective

2721:   Input Parameters:
2722: + options - options database, use `NULL` for default global database
2723: . pre     - string to prepend to each name or `NULL`
2724: - name    - the option one is seeking

2726:   Output Parameters:
2727: + dvalue - the double value to return
2728: - set    - `PETSC_TRUE` if found, `PETSC_FALSE` if not found

2730:   Level: beginner

2732:   Notes:
2733:   Accepts the special values `determine`, `decide` and `unlimited`.

2735:   Accepts the deprecated value `default`

2737:   If the user does not supply the option `dvalue` is NOT changed. Thus
2738:   you should ALWAYS initialize `dvalue` if you access it without first checking that the `set` flag is true.

2740: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2741:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2742:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2743:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2744:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2745:           `PetscOptionsFList()`, `PetscOptionsEList()`
2746: @*/
2747: PetscErrorCode PetscOptionsGetReal(PetscOptions options, const char pre[], const char name[], PetscReal *dvalue, PetscBool *set)
2748: {
2749:   const char *value;
2750:   PetscBool   flag;

2752:   PetscFunctionBegin;
2753:   PetscAssertPointer(name, 3);
2754:   PetscAssertPointer(dvalue, 4);
2755:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2756:   if (flag) {
2757:     if (!value) {
2758:       if (set) *set = PETSC_FALSE;
2759:     } else {
2760:       if (set) *set = PETSC_TRUE;
2761:       PetscCall(PetscOptionsStringToReal(value, dvalue));
2762:     }
2763:   } else {
2764:     if (set) *set = PETSC_FALSE;
2765:   }
2766:   PetscFunctionReturn(PETSC_SUCCESS);
2767: }

2769: /*@C
2770:   PetscOptionsGetScalar - Gets the scalar value for a particular
2771:   option in the database.

2773:   Not Collective

2775:   Input Parameters:
2776: + options - options database, use `NULL` for default global database
2777: . pre     - string to prepend to each name or `NULL`
2778: - name    - the option one is seeking

2780:   Output Parameters:
2781: + dvalue - the scalar value to return
2782: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2784:   Level: beginner

2786:   Example Usage:
2787:   A complex number 2+3i must be specified with NO spaces

2789:   Note:
2790:   If the user does not supply the option `dvalue` is NOT changed. Thus
2791:   you should ALWAYS initialize `dvalue` if you access it without first checking if the `set` flag is true.

2793: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2794:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2795:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2796:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2797:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2798:           `PetscOptionsFList()`, `PetscOptionsEList()`
2799: @*/
2800: PetscErrorCode PetscOptionsGetScalar(PetscOptions options, const char pre[], const char name[], PetscScalar *dvalue, PetscBool *set)
2801: {
2802:   const char *value;
2803:   PetscBool   flag;

2805:   PetscFunctionBegin;
2806:   PetscAssertPointer(name, 3);
2807:   PetscAssertPointer(dvalue, 4);
2808:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2809:   if (flag) {
2810:     if (!value) {
2811:       if (set) *set = PETSC_FALSE;
2812:     } else {
2813: #if !defined(PETSC_USE_COMPLEX)
2814:       PetscCall(PetscOptionsStringToReal(value, dvalue));
2815: #else
2816:       PetscCall(PetscOptionsStringToScalar(value, dvalue));
2817: #endif
2818:       if (set) *set = PETSC_TRUE;
2819:     }
2820:   } else { /* flag */
2821:     if (set) *set = PETSC_FALSE;
2822:   }
2823:   PetscFunctionReturn(PETSC_SUCCESS);
2824: }

2826: /*@C
2827:   PetscOptionsGetString - Gets the string value for a particular option in
2828:   the database.

2830:   Not Collective

2832:   Input Parameters:
2833: + options - options database, use `NULL` for default global database
2834: . pre     - string to prepend to name or `NULL`
2835: . name    - the option one is seeking
2836: - len     - maximum length of the string including null termination

2838:   Output Parameters:
2839: + string - location to copy string
2840: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2842:   Level: beginner

2844:   Note:
2845:   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`

2847:   If the user does not use the option then `string` is not changed. Thus
2848:   you should ALWAYS initialize `string` if you access it without first checking that the `set` flag is true.

2850:   Fortran Notes:
2851:   The Fortran interface is slightly different from the C/C++
2852:   interface (len is not used).  Sample usage in Fortran follows
2853: .vb
2854:       character *20    string
2855:       PetscErrorCode   ierr
2856:       PetscBool        set
2857:       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2858: .ve

2860: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2861:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2862:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2863:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2864:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2865:           `PetscOptionsFList()`, `PetscOptionsEList()`
2866: @*/
2867: PetscErrorCode PetscOptionsGetString(PetscOptions options, const char pre[], const char name[], char string[], size_t len, PetscBool *set)
2868: {
2869:   const char *value;
2870:   PetscBool   flag;

2872:   PetscFunctionBegin;
2873:   PetscAssertPointer(name, 3);
2874:   PetscAssertPointer(string, 4);
2875:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2876:   if (!flag) {
2877:     if (set) *set = PETSC_FALSE;
2878:   } else {
2879:     if (set) *set = PETSC_TRUE;
2880:     if (value) PetscCall(PetscStrncpy(string, value, len));
2881:     else PetscCall(PetscArrayzero(string, len));
2882:   }
2883:   PetscFunctionReturn(PETSC_SUCCESS);
2884: }

2886: /*@C
2887:   PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2888:   option in the database.  The values must be separated with commas with no intervening spaces.

2890:   Not Collective

2892:   Input Parameters:
2893: + options - options database, use `NULL` for default global database
2894: . pre     - string to prepend to each name or `NULL`
2895: - name    - the option one is seeking

2897:   Output Parameters:
2898: + dvalue - the Boolean values to return
2899: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2900: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2902:   Level: beginner

2904:   Note:
2905:   TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE`. FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE`

2907: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2908:           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2909:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2910:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2911:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2912:           `PetscOptionsFList()`, `PetscOptionsEList()`
2913: @*/
2914: PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options, const char pre[], const char name[], PetscBool dvalue[], PetscInt *nmax, PetscBool *set)
2915: {
2916:   const char *svalue;
2917:   char       *value;
2918:   PetscInt    n = 0;
2919:   PetscBool   flag;
2920:   PetscToken  token;

2922:   PetscFunctionBegin;
2923:   PetscAssertPointer(name, 3);
2924:   PetscAssertPointer(dvalue, 4);
2925:   PetscAssertPointer(nmax, 5);

2927:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2928:   if (!flag || !svalue) {
2929:     if (set) *set = PETSC_FALSE;
2930:     *nmax = 0;
2931:     PetscFunctionReturn(PETSC_SUCCESS);
2932:   }
2933:   if (set) *set = PETSC_TRUE;
2934:   PetscCall(PetscTokenCreate(svalue, ',', &token));
2935:   PetscCall(PetscTokenFind(token, &value));
2936:   while (value && n < *nmax) {
2937:     PetscCall(PetscOptionsStringToBool(value, dvalue));
2938:     PetscCall(PetscTokenFind(token, &value));
2939:     dvalue++;
2940:     n++;
2941:   }
2942:   PetscCall(PetscTokenDestroy(&token));
2943:   *nmax = n;
2944:   PetscFunctionReturn(PETSC_SUCCESS);
2945: }

2947: /*@C
2948:   PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.

2950:   Not Collective

2952:   Input Parameters:
2953: + options - options database, use `NULL` for default global database
2954: . pre     - option prefix or `NULL`
2955: . name    - option name
2956: - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2958:   Output Parameters:
2959: + ivalue - the  enum values to return
2960: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2961: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2963:   Level: beginner

2965:   Notes:
2966:   The array must be passed as a comma separated list with no spaces between the items.

2968:   `list` is usually something like `PCASMTypes` or some other predefined list of enum names.

2970: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2971:           `PetscOptionsGetEnum()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2972:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsName()`,
2973:           `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, `PetscOptionsStringArray()`, `PetscOptionsRealArray()`,
2974:           `PetscOptionsScalar()`, `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2975:           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2976: @*/
2977: PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options, const char pre[], const char name[], const char *const *list, PetscEnum ivalue[], PetscInt *nmax, PetscBool *set)
2978: {
2979:   const char *svalue;
2980:   char       *value;
2981:   PetscInt    n = 0;
2982:   PetscEnum   evalue;
2983:   PetscBool   flag;
2984:   PetscToken  token;

2986:   PetscFunctionBegin;
2987:   PetscAssertPointer(name, 3);
2988:   PetscAssertPointer(list, 4);
2989:   PetscAssertPointer(ivalue, 5);
2990:   PetscAssertPointer(nmax, 6);

2992:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2993:   if (!flag || !svalue) {
2994:     if (set) *set = PETSC_FALSE;
2995:     *nmax = 0;
2996:     PetscFunctionReturn(PETSC_SUCCESS);
2997:   }
2998:   if (set) *set = PETSC_TRUE;
2999:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3000:   PetscCall(PetscTokenFind(token, &value));
3001:   while (value && n < *nmax) {
3002:     PetscCall(PetscEnumFind(list, value, &evalue, &flag));
3003:     PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown enum value '%s' for -%s%s", svalue, pre ? pre : "", name + 1);
3004:     ivalue[n++] = evalue;
3005:     PetscCall(PetscTokenFind(token, &value));
3006:   }
3007:   PetscCall(PetscTokenDestroy(&token));
3008:   *nmax = n;
3009:   PetscFunctionReturn(PETSC_SUCCESS);
3010: }

3012: /*@C
3013:   PetscOptionsGetIntArray - Gets an array of integer values for a particular option in the database.

3015:   Not Collective

3017:   Input Parameters:
3018: + options - options database, use `NULL` for default global database
3019: . pre     - string to prepend to each name or `NULL`
3020: - name    - the option one is seeking

3022:   Output Parameters:
3023: + ivalue - the integer values to return
3024: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3025: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3027:   Level: beginner

3029:   Notes:
3030:   The array can be passed as
3031: +  a comma separated list -                                 0,1,2,3,4,5,6,7
3032: .  a range (start\-end+1) -                                 0-8
3033: .  a range with given increment (start\-end+1:inc) -        0-7:2
3034: -  a combination of values and ranges separated by commas - 0,1-8,8-15:2

3036:   There must be no intervening spaces between the values.

3038: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3039:           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3040:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3041:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3042:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3043:           `PetscOptionsFList()`, `PetscOptionsEList()`
3044: @*/
3045: PetscErrorCode PetscOptionsGetIntArray(PetscOptions options, const char pre[], const char name[], PetscInt ivalue[], PetscInt *nmax, PetscBool *set)
3046: {
3047:   const char *svalue;
3048:   char       *value;
3049:   PetscInt    n = 0, i, j, start, end, inc, nvalues;
3050:   size_t      len;
3051:   PetscBool   flag, foundrange;
3052:   PetscToken  token;

3054:   PetscFunctionBegin;
3055:   PetscAssertPointer(name, 3);
3056:   PetscAssertPointer(ivalue, 4);
3057:   PetscAssertPointer(nmax, 5);

3059:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3060:   if (!flag || !svalue) {
3061:     if (set) *set = PETSC_FALSE;
3062:     *nmax = 0;
3063:     PetscFunctionReturn(PETSC_SUCCESS);
3064:   }
3065:   if (set) *set = PETSC_TRUE;
3066:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3067:   PetscCall(PetscTokenFind(token, &value));
3068:   while (value && n < *nmax) {
3069:     /* look for form  d-D where d and D are integers */
3070:     foundrange = PETSC_FALSE;
3071:     PetscCall(PetscStrlen(value, &len));
3072:     if (value[0] == '-') i = 2;
3073:     else i = 1;
3074:     for (; i < (int)len; i++) {
3075:       if (value[i] == '-') {
3076:         PetscCheck(i != (int)len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value);
3077:         value[i] = 0;

3079:         PetscCall(PetscOptionsStringToInt(value, &start));
3080:         inc = 1;
3081:         j   = i + 1;
3082:         for (; j < (int)len; j++) {
3083:           if (value[j] == ':') {
3084:             value[j] = 0;

3086:             PetscCall(PetscOptionsStringToInt(value + j + 1, &inc));
3087:             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);
3088:             break;
3089:           }
3090:         }
3091:         PetscCall(PetscOptionsStringToInt(value + i + 1, &end));
3092:         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);
3093:         nvalues = (end - start) / inc + (end - start) % inc;
3094:         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);
3095:         for (; start < end; start += inc) {
3096:           *ivalue = start;
3097:           ivalue++;
3098:           n++;
3099:         }
3100:         foundrange = PETSC_TRUE;
3101:         break;
3102:       }
3103:     }
3104:     if (!foundrange) {
3105:       PetscCall(PetscOptionsStringToInt(value, ivalue));
3106:       ivalue++;
3107:       n++;
3108:     }
3109:     PetscCall(PetscTokenFind(token, &value));
3110:   }
3111:   PetscCall(PetscTokenDestroy(&token));
3112:   *nmax = n;
3113:   PetscFunctionReturn(PETSC_SUCCESS);
3114: }

3116: /*@C
3117:   PetscOptionsGetRealArray - Gets an array of double precision values for a
3118:   particular option in the database.  The values must be separated with commas with no intervening spaces.

3120:   Not Collective

3122:   Input Parameters:
3123: + options - options database, use `NULL` for default global database
3124: . pre     - string to prepend to each name or `NULL`
3125: - name    - the option one is seeking

3127:   Output Parameters:
3128: + dvalue - the double values to return
3129: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3130: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3132:   Level: beginner

3134: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3135:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3136:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3137:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3138:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3139:           `PetscOptionsFList()`, `PetscOptionsEList()`
3140: @*/
3141: PetscErrorCode PetscOptionsGetRealArray(PetscOptions options, const char pre[], const char name[], PetscReal dvalue[], PetscInt *nmax, PetscBool *set)
3142: {
3143:   const char *svalue;
3144:   char       *value;
3145:   PetscInt    n = 0;
3146:   PetscBool   flag;
3147:   PetscToken  token;

3149:   PetscFunctionBegin;
3150:   PetscAssertPointer(name, 3);
3151:   PetscAssertPointer(dvalue, 4);
3152:   PetscAssertPointer(nmax, 5);

3154:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3155:   if (!flag || !svalue) {
3156:     if (set) *set = PETSC_FALSE;
3157:     *nmax = 0;
3158:     PetscFunctionReturn(PETSC_SUCCESS);
3159:   }
3160:   if (set) *set = PETSC_TRUE;
3161:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3162:   PetscCall(PetscTokenFind(token, &value));
3163:   while (value && n < *nmax) {
3164:     PetscCall(PetscOptionsStringToReal(value, dvalue++));
3165:     PetscCall(PetscTokenFind(token, &value));
3166:     n++;
3167:   }
3168:   PetscCall(PetscTokenDestroy(&token));
3169:   *nmax = n;
3170:   PetscFunctionReturn(PETSC_SUCCESS);
3171: }

3173: /*@C
3174:   PetscOptionsGetScalarArray - Gets an array of scalars for a
3175:   particular option in the database.  The values must be separated with commas with no intervening spaces.

3177:   Not Collective

3179:   Input Parameters:
3180: + options - options database, use `NULL` for default global database
3181: . pre     - string to prepend to each name or `NULL`
3182: - name    - the option one is seeking

3184:   Output Parameters:
3185: + dvalue - the scalar values to return
3186: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3187: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3189:   Level: beginner

3191: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3192:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3193:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3194:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3195:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3196:           `PetscOptionsFList()`, `PetscOptionsEList()`
3197: @*/
3198: PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options, const char pre[], const char name[], PetscScalar dvalue[], PetscInt *nmax, PetscBool *set)
3199: {
3200:   const char *svalue;
3201:   char       *value;
3202:   PetscInt    n = 0;
3203:   PetscBool   flag;
3204:   PetscToken  token;

3206:   PetscFunctionBegin;
3207:   PetscAssertPointer(name, 3);
3208:   PetscAssertPointer(dvalue, 4);
3209:   PetscAssertPointer(nmax, 5);

3211:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3212:   if (!flag || !svalue) {
3213:     if (set) *set = PETSC_FALSE;
3214:     *nmax = 0;
3215:     PetscFunctionReturn(PETSC_SUCCESS);
3216:   }
3217:   if (set) *set = PETSC_TRUE;
3218:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3219:   PetscCall(PetscTokenFind(token, &value));
3220:   while (value && n < *nmax) {
3221:     PetscCall(PetscOptionsStringToScalar(value, dvalue++));
3222:     PetscCall(PetscTokenFind(token, &value));
3223:     n++;
3224:   }
3225:   PetscCall(PetscTokenDestroy(&token));
3226:   *nmax = n;
3227:   PetscFunctionReturn(PETSC_SUCCESS);
3228: }

3230: /*@C
3231:   PetscOptionsGetStringArray - Gets an array of string values for a particular
3232:   option in the database. The values must be separated with commas with no intervening spaces.

3234:   Not Collective; No Fortran Support

3236:   Input Parameters:
3237: + options - options database, use `NULL` for default global database
3238: . pre     - string to prepend to name or `NULL`
3239: - name    - the option one is seeking

3241:   Output Parameters:
3242: + strings - location to copy strings
3243: . nmax    - On input maximum number of strings, on output the actual number of strings found
3244: - set     - `PETSC_TRUE` if found, else `PETSC_FALSE`

3246:   Level: beginner

3248:   Notes:
3249:   The `nmax` parameter is used for both input and output.

3251:   The user should pass in an array of pointers to char, to hold all the
3252:   strings returned by this function.

3254:   The user is responsible for deallocating the strings that are
3255:   returned.

3257: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
3258:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3259:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3260:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3261:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3262:           `PetscOptionsFList()`, `PetscOptionsEList()`
3263: @*/
3264: PetscErrorCode PetscOptionsGetStringArray(PetscOptions options, const char pre[], const char name[], char *strings[], PetscInt *nmax, PetscBool *set)
3265: {
3266:   const char *svalue;
3267:   char       *value;
3268:   PetscInt    n = 0;
3269:   PetscBool   flag;
3270:   PetscToken  token;

3272:   PetscFunctionBegin;
3273:   PetscAssertPointer(name, 3);
3274:   PetscAssertPointer(strings, 4);
3275:   PetscAssertPointer(nmax, 5);

3277:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3278:   if (!flag || !svalue) {
3279:     if (set) *set = PETSC_FALSE;
3280:     *nmax = 0;
3281:     PetscFunctionReturn(PETSC_SUCCESS);
3282:   }
3283:   if (set) *set = PETSC_TRUE;
3284:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3285:   PetscCall(PetscTokenFind(token, &value));
3286:   while (value && n < *nmax) {
3287:     PetscCall(PetscStrallocpy(value, &strings[n]));
3288:     PetscCall(PetscTokenFind(token, &value));
3289:     n++;
3290:   }
3291:   PetscCall(PetscTokenDestroy(&token));
3292:   *nmax = n;
3293:   PetscFunctionReturn(PETSC_SUCCESS);
3294: }

3296: /*@C
3297:   PetscOptionsDeprecated_Private - mark an option as deprecated, optionally replacing it with `newname`

3299:   Prints a deprecation warning, unless an option is supplied to suppress.

3301:   Logically Collective

3303:   Input Parameters:
3304: + PetscOptionsObject - string to prepend to name or `NULL`
3305: . oldname            - the old, deprecated option
3306: . newname            - the new option, or `NULL` if option is purely removed
3307: . version            - a string describing the version of first deprecation, e.g. "3.9"
3308: - info               - additional information string, or `NULL`.

3310:   Options Database Key:
3311: . -options_suppress_deprecated_warnings - do not print deprecation warnings

3313:   Level: developer

3315:   Notes:
3316:   If `newname` is provided then the options database will automatically check the database for `oldname`.

3318:   The old call `PetscOptionsXXX`(`oldname`) should be removed from the source code when both (1) the call to `PetscOptionsDeprecated()` occurs before the
3319:   new call to `PetscOptionsXXX`(`newname`) and (2) the argument handling of the new call to `PetscOptionsXXX`(`newname`) is identical to the previous call.
3320:   See `PTScotch_PartGraph_Seq()` for an example of when (1) fails and `SNESTestJacobian()` where an example of (2) fails.

3322:   Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
3323:   Only the process of rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or
3324:   `PetscObjectOptionsBegin()` prints the information
3325:   If newname is provided, the old option is replaced. Otherwise, it remains
3326:   in the options database.
3327:   If an option is not replaced, the info argument should be used to advise the user
3328:   on how to proceed.
3329:   There is a limit on the length of the warning printed, so very long strings
3330:   provided as info may be truncated.

3332: .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
3333: @*/
3334: PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject, const char oldname[], const char newname[], const char version[], const char info[])
3335: {
3336:   PetscBool         found, quiet;
3337:   const char       *value;
3338:   const char *const quietopt = "-options_suppress_deprecated_warnings";
3339:   char              msg[4096];
3340:   char             *prefix  = NULL;
3341:   PetscOptions      options = NULL;
3342:   MPI_Comm          comm    = PETSC_COMM_SELF;

3344:   PetscFunctionBegin;
3345:   PetscAssertPointer(oldname, 2);
3346:   PetscAssertPointer(version, 4);
3347:   if (PetscOptionsObject) {
3348:     prefix  = PetscOptionsObject->prefix;
3349:     options = PetscOptionsObject->options;
3350:     comm    = PetscOptionsObject->comm;
3351:   }
3352:   PetscCall(PetscOptionsFindPair(options, prefix, oldname, &value, &found));
3353:   if (found) {
3354:     if (newname) {
3355:       if (prefix) PetscCall(PetscOptionsPrefixPush(options, prefix));
3356:       PetscCall(PetscOptionsSetValue(options, newname, value));
3357:       if (prefix) PetscCall(PetscOptionsPrefixPop(options));
3358:       PetscCall(PetscOptionsClearValue(options, oldname));
3359:     }
3360:     quiet = PETSC_FALSE;
3361:     PetscCall(PetscOptionsGetBool(options, NULL, quietopt, &quiet, NULL));
3362:     if (!quiet) {
3363:       PetscCall(PetscStrncpy(msg, "** PETSc DEPRECATION WARNING ** : the option ", sizeof(msg)));
3364:       PetscCall(PetscStrlcat(msg, oldname, sizeof(msg)));
3365:       PetscCall(PetscStrlcat(msg, " is deprecated as of version ", sizeof(msg)));
3366:       PetscCall(PetscStrlcat(msg, version, sizeof(msg)));
3367:       PetscCall(PetscStrlcat(msg, " and will be removed in a future release.\n", sizeof(msg)));
3368:       if (newname) {
3369:         PetscCall(PetscStrlcat(msg, "   Use the option ", sizeof(msg)));
3370:         PetscCall(PetscStrlcat(msg, newname, sizeof(msg)));
3371:         PetscCall(PetscStrlcat(msg, " instead.", sizeof(msg)));
3372:       }
3373:       if (info) {
3374:         PetscCall(PetscStrlcat(msg, " ", sizeof(msg)));
3375:         PetscCall(PetscStrlcat(msg, info, sizeof(msg)));
3376:       }
3377:       PetscCall(PetscStrlcat(msg, " (Silence this warning with ", sizeof(msg)));
3378:       PetscCall(PetscStrlcat(msg, quietopt, sizeof(msg)));
3379:       PetscCall(PetscStrlcat(msg, ")\n", sizeof(msg)));
3380:       PetscCall(PetscPrintf(comm, "%s", msg));
3381:     }
3382:   }
3383:   PetscFunctionReturn(PETSC_SUCCESS);
3384: }