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:   PetscCtxDestroyFn *monitordestroy[MAXOPTIONSMONITORS];                                                /* 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, fname, sizeof(fname)));
473:     PetscCall(PetscFixFilename(fname, fpath));
474:     PetscCall(PetscGetFullPath(fpath, fname, sizeof(fname)));

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

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

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

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

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

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

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

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

602:   Collective

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

612:   Level: developer

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

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

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

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

649:   Logically Collective

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

656:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

809:   Collective on `PETSC_COMM_WORLD`

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

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

824:   Level: advanced

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

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

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

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

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

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

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

888:   /* insert YAML environment options */
889:   if (rank == 0) {
890:     eoptions = getenv("PETSC_OPTIONS_YAML");
891:     PetscCall(PetscStrlen(eoptions, &len));
892:   }
893:   PetscCallMPI(MPI_Bcast(&len, 1, MPIU_SIZE_T, 0, comm));
894:   if (len) {
895:     if (rank) PetscCall(PetscMalloc1(len + 1, &eoptions));
896:     PetscCallMPI(MPI_Bcast(eoptions, (PetscMPIInt)len, MPI_CHAR, 0, comm));
897:     if (rank) eoptions[len] = 0;
898:     PetscCall(PetscOptionsInsertStringYAML_Private(options, eoptions, PETSC_OPT_ENVIRONMENT));
899:     if (rank) PetscCall(PetscFree(eoptions));
900:   }

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

908: /* These options are not printed with PetscOptionsView() or PetscOptionsMonitor() when PetscCIEnabled is on */
909: /* TODO: get the list from the test harness, do not have it hardwired here. Maybe from gmakegentest.py */
910: static const char *PetscCIOptions[] = {"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"};

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

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

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

926:   Logically Collective, No Fortran Support

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

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

935:   Level: advanced

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

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

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

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

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

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

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

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

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

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

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

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

1034:   Logically Collective

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

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

1044:   Level: advanced

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

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

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

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

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

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

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

1092:   Level: advanced

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

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

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

1112:   Logically Collective

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

1117:   Level: developer

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

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

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

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

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

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

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

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

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

1177:   Logically Collective

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

1184:   Level: advanced

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

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

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

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

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

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

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

1247:   Logically Collective

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

1254:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1405:   Logically Collective

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

1411:   Level: intermediate

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

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

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

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

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

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

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

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

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

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

1474:   Not Collective

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

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

1485:   Level: developer

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

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

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

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

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

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

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

1543:   khash_t(HO) *ht = options->ht;
1544:   khiter_t     it = kh_get(HO, ht, name);
1545:   if (it != kh_end(ht)) {
1546:     int i            = kh_val(ht, it);
1547:     options->used[i] = PETSC_TRUE;
1548:     if (value) *value = options->values[i];
1549:     if (set) *set = PETSC_TRUE;
1550:     PetscFunctionReturn(PETSC_SUCCESS);
1551:   }

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

1588:   if (set) *set = PETSC_FALSE;
1589:   PetscFunctionReturn(PETSC_SUCCESS);
1590: }

1592: /* Check whether any option begins with pre+name */
1593: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options, const char pre[], const char name[], const char *option[], const char *value[], PetscBool *set)
1594: {
1595:   char buf[PETSC_MAX_OPTION_NAME];
1596:   int  numCnt = 0, locs[16], loce[16];

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

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

1605:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1606:   if (pre && pre[0]) {
1607:     char *ptr = buf;
1608:     if (name[0] == '-') {
1609:       *ptr++ = '-';
1610:       name++;
1611:     }
1612:     PetscCall(PetscStrncpy(ptr, pre, sizeof(buf) - ((ptr == buf) ? 0 : 1)));
1613:     PetscCall(PetscStrlcat(buf, name, sizeof(buf)));
1614:     name = buf;
1615:   }

1617:   if (PetscDefined(USE_DEBUG)) {
1618:     PetscBool valid;
1619:     char      key[PETSC_MAX_OPTION_NAME + 1] = "-";
1620:     PetscCall(PetscStrncpy(key + 1, name, sizeof(key) - 1));
1621:     PetscCall(PetscOptionsValidKey(key, &valid));
1622:     PetscCheck(valid, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid option '%s' obtained from pre='%s' and name='%s'", key, pre ? pre : "", name);
1623:   }

1625:   /* determine the location and number of all _%d_ in the key */
1626:   {
1627:     int i, j;
1628:     for (i = 0; name[i]; i++) {
1629:       if (name[i] == '_') {
1630:         for (j = i + 1; name[j]; j++) {
1631:           if (name[j] >= '0' && name[j] <= '9') continue;
1632:           if (name[j] == '_' && j > i + 1) { /* found a number */
1633:             locs[numCnt]   = i + 1;
1634:             loce[numCnt++] = j + 1;
1635:           }
1636:           i = j - 1;
1637:           break;
1638:         }
1639:       }
1640:     }
1641:   }

1643:   /* slow search */
1644:   for (int c = -1; c < numCnt; ++c) {
1645:     char   opt[PETSC_MAX_OPTION_NAME + 2] = "";
1646:     size_t len;

1648:     if (c < 0) {
1649:       PetscCall(PetscStrncpy(opt, name, sizeof(opt)));
1650:     } else {
1651:       PetscCall(PetscStrncpy(opt, name, PetscMin((size_t)(locs[c] + 1), sizeof(opt))));
1652:       PetscCall(PetscStrlcat(opt, name + loce[c], sizeof(opt) - 1));
1653:     }
1654:     PetscCall(PetscStrlen(opt, &len));
1655:     for (int i = 0; i < options->N; i++) {
1656:       PetscBool match;

1658:       PetscCall(PetscStrncmp(options->names[i], opt, len, &match));
1659:       if (match) {
1660:         options->used[i] = PETSC_TRUE;
1661:         if (option) *option = options->names[i];
1662:         if (value) *value = options->values[i];
1663:         if (set) *set = PETSC_TRUE;
1664:         PetscFunctionReturn(PETSC_SUCCESS);
1665:       }
1666:     }
1667:   }

1669:   if (set) *set = PETSC_FALSE;
1670:   PetscFunctionReturn(PETSC_SUCCESS);
1671: }

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

1676:   Not Collective

1678:   Input Parameters:
1679: + options - options database, use `NULL` for default global database
1680: . pre     - the option prefix (may be `NULL`)
1681: . name    - the option name one is seeking
1682: - mess    - error message (may be `NULL`)

1684:   Level: advanced

1686: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`,
1687:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1688:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1689:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1690:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1691:           `PetscOptionsFList()`, `PetscOptionsEList()`
1692: @*/
1693: PetscErrorCode PetscOptionsReject(PetscOptions options, const char pre[], const char name[], const char mess[])
1694: {
1695:   PetscBool flag = PETSC_FALSE;

1697:   PetscFunctionBegin;
1698:   PetscCall(PetscOptionsHasName(options, pre, name, &flag));
1699:   if (flag) {
1700:     PetscCheck(!mess || !mess[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s with %s", pre ? pre : "", name + 1, mess);
1701:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Program has disabled option: -%s%s", pre ? pre : "", name + 1);
1702:   }
1703:   PetscFunctionReturn(PETSC_SUCCESS);
1704: }

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

1709:   Not Collective

1711:   Input Parameter:
1712: . options - options database, use `NULL` for default global database

1714:   Output Parameter:
1715: . set - `PETSC_TRUE` if found else `PETSC_FALSE`.

1717:   Level: advanced

1719: .seealso: `PetscOptionsHasName()`
1720: @*/
1721: PetscErrorCode PetscOptionsHasHelp(PetscOptions options, PetscBool *set)
1722: {
1723:   PetscFunctionBegin;
1724:   PetscAssertPointer(set, 2);
1725:   options = options ? options : defaultoptions;
1726:   *set    = options->help;
1727:   PetscFunctionReturn(PETSC_SUCCESS);
1728: }

1730: PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options, PetscBool *set)
1731: {
1732:   PetscFunctionBegin;
1733:   PetscAssertPointer(set, 2);
1734:   options = options ? options : defaultoptions;
1735:   *set    = options->help_intro;
1736:   PetscFunctionReturn(PETSC_SUCCESS);
1737: }

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

1743:   Not Collective

1745:   Input Parameters:
1746: + options - options database, use `NULL` for default global database
1747: . pre     - string to prepend to the name or `NULL`
1748: - name    - the option one is seeking

1750:   Output Parameter:
1751: . set - `PETSC_TRUE` if found else `PETSC_FALSE`.

1753:   Level: beginner

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

1758: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
1759:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
1760:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1761:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1762:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1763:           `PetscOptionsFList()`, `PetscOptionsEList()`
1764: @*/
1765: PetscErrorCode PetscOptionsHasName(PetscOptions options, const char pre[], const char name[], PetscBool *set)
1766: {
1767:   const char *value;
1768:   PetscBool   flag;

1770:   PetscFunctionBegin;
1771:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
1772:   if (set) *set = flag;
1773:   PetscFunctionReturn(PETSC_SUCCESS);
1774: }

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

1779:   Not Collective

1781:   Input Parameter:
1782: . options - the options database, use `NULL` for the default global database

1784:   Output Parameter:
1785: . copts - pointer where string pointer is stored

1787:   Level: advanced

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

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

1794: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsView()`, `PetscOptionsPush()`, `PetscOptionsPop()`
1795: @*/
1796: PetscErrorCode PetscOptionsGetAll(PetscOptions options, char *copts[])
1797: {
1798:   PetscInt i;
1799:   size_t   len = 1, lent = 0;
1800:   char    *coptions = NULL;

1802:   PetscFunctionBegin;
1803:   PetscAssertPointer(copts, 2);
1804:   options = options ? options : defaultoptions;
1805:   /* count the length of the required string */
1806:   for (i = 0; i < options->N; i++) {
1807:     PetscCall(PetscStrlen(options->names[i], &lent));
1808:     len += 2 + lent;
1809:     if (options->values[i]) {
1810:       PetscCall(PetscStrlen(options->values[i], &lent));
1811:       len += 1 + lent;
1812:     }
1813:   }
1814:   PetscCall(PetscMalloc1(len, &coptions));
1815:   coptions[0] = 0;
1816:   for (i = 0; i < options->N; i++) {
1817:     PetscCall(PetscStrlcat(coptions, "-", len));
1818:     PetscCall(PetscStrlcat(coptions, options->names[i], len));
1819:     PetscCall(PetscStrlcat(coptions, " ", len));
1820:     if (options->values[i]) {
1821:       PetscCall(PetscStrlcat(coptions, options->values[i], len));
1822:       PetscCall(PetscStrlcat(coptions, " ", len));
1823:     }
1824:   }
1825:   *copts = coptions;
1826:   PetscFunctionReturn(PETSC_SUCCESS);
1827: }

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

1832:   Not Collective

1834:   Input Parameters:
1835: + options - options database, use `NULL` for default global database
1836: - name    - string name of option

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

1841:   Level: advanced

1843:   Note:
1844:   The value returned may be different on each process and depends on which options have been processed
1845:   on the given process

1847: .seealso: `PetscOptionsView()`, `PetscOptionsLeft()`, `PetscOptionsAllUsed()`
1848: @*/
1849: PetscErrorCode PetscOptionsUsed(PetscOptions options, const char *name, PetscBool *used)
1850: {
1851:   PetscInt i;

1853:   PetscFunctionBegin;
1854:   PetscAssertPointer(name, 2);
1855:   PetscAssertPointer(used, 3);
1856:   options = options ? options : defaultoptions;
1857:   *used   = PETSC_FALSE;
1858:   for (i = 0; i < options->N; i++) {
1859:     PetscCall(PetscStrcasecmp(options->names[i], name, used));
1860:     if (*used) {
1861:       *used = options->used[i];
1862:       break;
1863:     }
1864:   }
1865:   PetscFunctionReturn(PETSC_SUCCESS);
1866: }

1868: /*@
1869:   PetscOptionsAllUsed - Returns a count of the number of options in the
1870:   database that have never been selected.

1872:   Not Collective

1874:   Input Parameter:
1875: . options - options database, use `NULL` for default global database

1877:   Output Parameter:
1878: . N - count of options not used

1880:   Level: advanced

1882:   Note:
1883:   The value returned may be different on each process and depends on which options have been processed
1884:   on the given process

1886: .seealso: `PetscOptionsView()`
1887: @*/
1888: PetscErrorCode PetscOptionsAllUsed(PetscOptions options, PetscInt *N)
1889: {
1890:   PetscInt i, n = 0;

1892:   PetscFunctionBegin;
1893:   PetscAssertPointer(N, 2);
1894:   options = options ? options : defaultoptions;
1895:   for (i = 0; i < options->N; i++) {
1896:     if (!options->used[i]) n++;
1897:   }
1898:   *N = n;
1899:   PetscFunctionReturn(PETSC_SUCCESS);
1900: }

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

1905:   Not Collective

1907:   Input Parameter:
1908: . options - options database; use `NULL` for default global database

1910:   Options Database Key:
1911: . -options_left - activates `PetscOptionsAllUsed()` within `PetscFinalize()`

1913:   Level: advanced

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

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

1923: .seealso: `PetscOptionsAllUsed()`
1924: @*/
1925: PetscErrorCode PetscOptionsLeft(PetscOptions options)
1926: {
1927:   PetscInt     i;
1928:   PetscInt     cnt = 0;
1929:   PetscOptions toptions;

1931:   PetscFunctionBegin;
1932:   toptions = options ? options : defaultoptions;
1933:   for (i = 0; i < toptions->N; i++) {
1934:     if (!toptions->used[i]) {
1935:       if (PetscCIOption(toptions->names[i])) continue;
1936:       if (toptions->values[i]) {
1937:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s value: %s source: %s\n", toptions->names[i], toptions->values[i], PetscOptionSources[toptions->source[i]]));
1938:       } else {
1939:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Option left: name:-%s (no value) source: %s\n", toptions->names[i], PetscOptionSources[toptions->source[i]]));
1940:       }
1941:     }
1942:   }
1943:   if (!options) {
1944:     toptions = defaultoptions;
1945:     while (toptions->previous) {
1946:       cnt++;
1947:       toptions = toptions->previous;
1948:     }
1949:     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));
1950:   }
1951:   PetscFunctionReturn(PETSC_SUCCESS);
1952: }

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

1957:   Not Collective

1959:   Input Parameter:
1960: . options - options database, use `NULL` for default global database

1962:   Output Parameters:
1963: + N      - count of options not used
1964: . names  - names of options not used
1965: - values - values of options not used

1967:   Level: advanced

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

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

1975: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`
1976: @*/
1977: PetscErrorCode PetscOptionsLeftGet(PetscOptions options, PetscInt *N, char **names[], char **values[])
1978: {
1979:   PetscInt i, n;

1981:   PetscFunctionBegin;
1982:   if (N) PetscAssertPointer(N, 2);
1983:   if (names) PetscAssertPointer(names, 3);
1984:   if (values) PetscAssertPointer(values, 4);
1985:   options = options ? options : defaultoptions;

1987:   /* The number of unused PETSc options */
1988:   n = 0;
1989:   for (i = 0; i < options->N; i++) {
1990:     if (PetscCIOption(options->names[i])) continue;
1991:     if (!options->used[i]) n++;
1992:   }
1993:   if (N) *N = n;
1994:   if (names) PetscCall(PetscMalloc1(n, names));
1995:   if (values) PetscCall(PetscMalloc1(n, values));

1997:   n = 0;
1998:   if (names || values) {
1999:     for (i = 0; i < options->N; i++) {
2000:       if (!options->used[i]) {
2001:         if (PetscCIOption(options->names[i])) continue;
2002:         if (names) (*names)[n] = options->names[i];
2003:         if (values) (*values)[n] = options->values[i];
2004:         n++;
2005:       }
2006:     }
2007:   }
2008:   PetscFunctionReturn(PETSC_SUCCESS);
2009: }

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

2014:   Not Collective

2016:   Input Parameters:
2017: + options - options database, use `NULL` for default global database
2018: . N       - count of options not used
2019: . names   - names of options not used
2020: - values  - values of options not used

2022:   Level: advanced

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

2027: .seealso: `PetscOptionsAllUsed()`, `PetscOptionsLeft()`, `PetscOptionsLeftGet()`
2028: @*/
2029: PetscErrorCode PetscOptionsLeftRestore(PetscOptions options, PetscInt *N, char **names[], char **values[])
2030: {
2031:   PetscFunctionBegin;
2032:   (void)options;
2033:   if (N) PetscAssertPointer(N, 2);
2034:   if (names) PetscAssertPointer(names, 3);
2035:   if (values) PetscAssertPointer(values, 4);
2036:   if (N) *N = 0;
2037:   if (names) PetscCall(PetscFree(*names));
2038:   if (values) PetscCall(PetscFree(*values));
2039:   PetscFunctionReturn(PETSC_SUCCESS);
2040: }

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

2045:   Logically Collective

2047:   Input Parameters:
2048: + name   - option name string
2049: . value  - option value string
2050: . source - The source for the option
2051: - ctx    - a `PETSCVIEWERASCII` or `NULL`

2053:   Level: intermediate

2055:   Notes:
2056:   If ctx is `NULL`, `PetscPrintf()` is used.
2057:   The first MPI process in the `PetscViewer` viewer actually prints the values, other
2058:   processes may have different values set

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

2062: .seealso: `PetscOptionsMonitorSet()`
2063: @*/
2064: PetscErrorCode PetscOptionsMonitorDefault(const char name[], const char value[], PetscOptionSource source, void *ctx)
2065: {
2066:   PetscFunctionBegin;
2067:   if (PetscCIOption(name)) PetscFunctionReturn(PETSC_SUCCESS);

2069:   if (ctx) {
2070:     PetscViewer viewer = (PetscViewer)ctx;
2071:     if (!value) {
2072:       PetscCall(PetscViewerASCIIPrintf(viewer, "Removing option: %s\n", name));
2073:     } else if (!value[0]) {
2074:       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
2075:     } else {
2076:       PetscCall(PetscViewerASCIIPrintf(viewer, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
2077:     }
2078:   } else {
2079:     MPI_Comm comm = PETSC_COMM_WORLD;
2080:     if (!value) {
2081:       PetscCall(PetscPrintf(comm, "Removing option: %s\n", name));
2082:     } else if (!value[0]) {
2083:       PetscCall(PetscPrintf(comm, "Setting option: %s (no value) (source: %s)\n", name, PetscOptionSources[source]));
2084:     } else {
2085:       PetscCall(PetscPrintf(comm, "Setting option: %s = %s (source: %s)\n", name, value, PetscOptionSources[source]));
2086:     }
2087:   }
2088:   PetscFunctionReturn(PETSC_SUCCESS);
2089: }

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

2095:   Not Collective

2097:   Input Parameters:
2098: + monitor        - pointer to function (if this is `NULL`, it turns off monitoring
2099: . mctx           - [optional] context for private data for the monitor routine (use `NULL` if
2100:                    no context is desired)
2101: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for its calling sequence

2103:   Calling sequence of `monitor`:
2104: + name   - option name string
2105: . value  - option value string, a value of `NULL` indicates the option is being removed from the database. A value
2106:            of "" indicates the option is in the database but has no value.
2107: . source - option source
2108: - mctx   - optional monitoring context, as set by `PetscOptionsMonitorSet()`

2110:   Options Database Keys:
2111: + -options_monitor <viewer> - turn on default monitoring of changes to the options database
2112: - -options_monitor_cancel   - turn off any option monitors except the default monitor obtained with `-options_monitor`

2114:   Level: intermediate

2116:   Notes:
2117:   See `PetscInitialize()` for options related to option database monitoring.

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

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

2127: .seealso: `PetscOptionsMonitorDefault()`, `PetscInitialize()`, `PetscCtxDestroyFn`
2128: @*/
2129: PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], PetscOptionSource source, void *mctx), void *mctx, PetscCtxDestroyFn *monitordestroy)
2130: {
2131:   PetscOptions options = defaultoptions;

2133:   PetscFunctionBegin;
2134:   if (options->monitorCancel) PetscFunctionReturn(PETSC_SUCCESS);
2135:   PetscCheck(options->numbermonitors < MAXOPTIONSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many PetscOptions monitors set");
2136:   options->monitor[options->numbermonitors]          = monitor;
2137:   options->monitordestroy[options->numbermonitors]   = monitordestroy;
2138:   options->monitorcontext[options->numbermonitors++] = mctx;
2139:   PetscFunctionReturn(PETSC_SUCCESS);
2140: }

2142: /*
2143:    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
2144:      Empty string is considered as true.
2145: */
2146: PetscErrorCode PetscOptionsStringToBool(const char value[], PetscBool *a)
2147: {
2148:   PetscBool istrue, isfalse;
2149:   size_t    len;

2151:   PetscFunctionBegin;
2152:   /* PetscStrlen() returns 0 for NULL or "" */
2153:   PetscCall(PetscStrlen(value, &len));
2154:   if (!len) {
2155:     *a = PETSC_TRUE;
2156:     PetscFunctionReturn(PETSC_SUCCESS);
2157:   }
2158:   PetscCall(PetscStrcasecmp(value, "TRUE", &istrue));
2159:   if (istrue) {
2160:     *a = PETSC_TRUE;
2161:     PetscFunctionReturn(PETSC_SUCCESS);
2162:   }
2163:   PetscCall(PetscStrcasecmp(value, "YES", &istrue));
2164:   if (istrue) {
2165:     *a = PETSC_TRUE;
2166:     PetscFunctionReturn(PETSC_SUCCESS);
2167:   }
2168:   PetscCall(PetscStrcasecmp(value, "1", &istrue));
2169:   if (istrue) {
2170:     *a = PETSC_TRUE;
2171:     PetscFunctionReturn(PETSC_SUCCESS);
2172:   }
2173:   PetscCall(PetscStrcasecmp(value, "on", &istrue));
2174:   if (istrue) {
2175:     *a = PETSC_TRUE;
2176:     PetscFunctionReturn(PETSC_SUCCESS);
2177:   }
2178:   PetscCall(PetscStrcasecmp(value, "FALSE", &isfalse));
2179:   if (isfalse) {
2180:     *a = PETSC_FALSE;
2181:     PetscFunctionReturn(PETSC_SUCCESS);
2182:   }
2183:   PetscCall(PetscStrcasecmp(value, "NO", &isfalse));
2184:   if (isfalse) {
2185:     *a = PETSC_FALSE;
2186:     PetscFunctionReturn(PETSC_SUCCESS);
2187:   }
2188:   PetscCall(PetscStrcasecmp(value, "0", &isfalse));
2189:   if (isfalse) {
2190:     *a = PETSC_FALSE;
2191:     PetscFunctionReturn(PETSC_SUCCESS);
2192:   }
2193:   PetscCall(PetscStrcasecmp(value, "off", &isfalse));
2194:   if (isfalse) {
2195:     *a = PETSC_FALSE;
2196:     PetscFunctionReturn(PETSC_SUCCESS);
2197:   }
2198:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown logical value: %s", value);
2199: }

2201: /*
2202:    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2203: */
2204: PetscErrorCode PetscOptionsStringToInt(const char name[], PetscInt *a)
2205: {
2206:   size_t    len;
2207:   PetscBool decide, tdefault, mouse, unlimited;

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

2213:   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &tdefault));
2214:   if (!tdefault) PetscCall(PetscStrcasecmp(name, "DEFAULT", &tdefault));
2215:   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &decide));
2216:   if (!decide) PetscCall(PetscStrcasecmp(name, "DECIDE", &decide));
2217:   if (!decide) PetscCall(PetscStrcasecmp(name, "PETSC_DETERMINE", &decide));
2218:   if (!decide) PetscCall(PetscStrcasecmp(name, "DETERMINE", &decide));
2219:   PetscCall(PetscStrcasecmp(name, "PETSC_UNLIMITED", &unlimited));
2220:   if (!unlimited) PetscCall(PetscStrcasecmp(name, "UNLIMITED", &unlimited));
2221:   PetscCall(PetscStrcasecmp(name, "mouse", &mouse));

2223:   if (tdefault) *a = PETSC_DEFAULT;
2224:   else if (decide) *a = PETSC_DECIDE;
2225:   else if (unlimited) *a = PETSC_UNLIMITED;
2226:   else if (mouse) *a = -1;
2227:   else {
2228:     char *endptr;
2229:     long  strtolval;

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

2234: #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2235:     (void)strtolval;
2236:     *a = atoll(name);
2237: #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2238:     (void)strtolval;
2239:     *a = _atoi64(name);
2240: #else
2241:     *a = (PetscInt)strtolval;
2242: #endif
2243:   }
2244:   PetscFunctionReturn(PETSC_SUCCESS);
2245: }

2247: #if defined(PETSC_USE_REAL___FLOAT128)
2248:   #include <quadmath.h>
2249: #endif

2251: static PetscErrorCode PetscStrtod(const char name[], PetscReal *a, char **endptr)
2252: {
2253:   PetscFunctionBegin;
2254: #if defined(PETSC_USE_REAL___FLOAT128)
2255:   *a = strtoflt128(name, endptr);
2256: #else
2257:   *a = (PetscReal)strtod(name, endptr);
2258: #endif
2259:   PetscFunctionReturn(PETSC_SUCCESS);
2260: }

2262: static PetscErrorCode PetscStrtoz(const char name[], PetscScalar *a, char **endptr, PetscBool *isImaginary)
2263: {
2264:   PetscBool hasi = PETSC_FALSE;
2265:   char     *ptr;
2266:   PetscReal strtoval;

2268:   PetscFunctionBegin;
2269:   PetscCall(PetscStrtod(name, &strtoval, &ptr));
2270:   if (ptr == name) {
2271:     strtoval = 1.;
2272:     hasi     = PETSC_TRUE;
2273:     if (name[0] == 'i') {
2274:       ptr++;
2275:     } else if (name[0] == '+' && name[1] == 'i') {
2276:       ptr += 2;
2277:     } else if (name[0] == '-' && name[1] == 'i') {
2278:       strtoval = -1.;
2279:       ptr += 2;
2280:     }
2281:   } else if (*ptr == 'i') {
2282:     hasi = PETSC_TRUE;
2283:     ptr++;
2284:   }
2285:   *endptr      = ptr;
2286:   *isImaginary = hasi;
2287:   if (hasi) {
2288: #if !defined(PETSC_USE_COMPLEX)
2289:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s contains imaginary but complex not supported ", name);
2290: #else
2291:     *a = PetscCMPLX(0., strtoval);
2292: #endif
2293:   } else {
2294:     *a = strtoval;
2295:   }
2296:   PetscFunctionReturn(PETSC_SUCCESS);
2297: }

2299: /*
2300:    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2301: */
2302: PetscErrorCode PetscOptionsStringToReal(const char name[], PetscReal *a)
2303: {
2304:   size_t    len;
2305:   PetscBool match;
2306:   char     *endptr;

2308:   PetscFunctionBegin;
2309:   PetscCall(PetscStrlen(name, &len));
2310:   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "String of length zero has no numerical value");

2312:   PetscCall(PetscStrcasecmp(name, "PETSC_DEFAULT", &match));
2313:   if (!match) PetscCall(PetscStrcasecmp(name, "DEFAULT", &match));
2314:   if (match) {
2315:     *a = PETSC_DEFAULT;
2316:     PetscFunctionReturn(PETSC_SUCCESS);
2317:   }

2319:   PetscCall(PetscStrcasecmp(name, "PETSC_DECIDE", &match));
2320:   if (!match) PetscCall(PetscStrcasecmp(name, "DECIDE", &match));
2321:   if (match) {
2322:     *a = PETSC_DECIDE;
2323:     PetscFunctionReturn(PETSC_SUCCESS);
2324:   }

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

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

2340:   PetscCall(PetscStrtod(name, a, &endptr));
2341:   PetscCheck((size_t)(endptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value", name);
2342:   PetscFunctionReturn(PETSC_SUCCESS);
2343: }

2345: PetscErrorCode PetscOptionsStringToScalar(const char name[], PetscScalar *a)
2346: {
2347:   PetscBool   imag1;
2348:   size_t      len;
2349:   PetscScalar val = 0.;
2350:   char       *ptr = NULL;

2352:   PetscFunctionBegin;
2353:   PetscCall(PetscStrlen(name, &len));
2354:   PetscCheck(len, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "character string of length zero has no numerical value");
2355:   PetscCall(PetscStrtoz(name, &val, &ptr, &imag1));
2356: #if defined(PETSC_USE_COMPLEX)
2357:   if ((size_t)(ptr - name) < len) {
2358:     PetscBool   imag2;
2359:     PetscScalar val2;

2361:     PetscCall(PetscStrtoz(ptr, &val2, &ptr, &imag2));
2362:     if (imag1) PetscCheck(imag2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s: must specify imaginary component second", name);
2363:     val = PetscCMPLX(PetscRealPart(val), PetscImaginaryPart(val2));
2364:   }
2365: #endif
2366:   PetscCheck((size_t)(ptr - name) == len, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input string %s has no numeric value ", name);
2367:   *a = val;
2368:   PetscFunctionReturn(PETSC_SUCCESS);
2369: }

2371: /*@C
2372:   PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2373:   option in the database.

2375:   Not Collective

2377:   Input Parameters:
2378: + options - options database, use `NULL` for default global database
2379: . pre     - the string to prepend to the name or `NULL`
2380: - name    - the option one is seeking

2382:   Output Parameters:
2383: + ivalue - the logical value to return
2384: - set    - `PETSC_TRUE`  if found, else `PETSC_FALSE`

2386:   Level: beginner

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

2392:   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`
2393:   is equivalent to `-requested_bool true`

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

2398: .seealso: `PetscOptionsGetBool3()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2399:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2400:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2401:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2402:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2403:           `PetscOptionsFList()`, `PetscOptionsEList()`
2404: @*/
2405: PetscErrorCode PetscOptionsGetBool(PetscOptions options, const char pre[], const char name[], PetscBool *ivalue, PetscBool *set)
2406: {
2407:   const char *value;
2408:   PetscBool   flag;

2410:   PetscFunctionBegin;
2411:   PetscAssertPointer(name, 3);
2412:   if (ivalue) PetscAssertPointer(ivalue, 4);
2413:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2414:   if (flag) {
2415:     if (set) *set = PETSC_TRUE;
2416:     PetscCall(PetscOptionsStringToBool(value, &flag));
2417:     if (ivalue) *ivalue = flag;
2418:   } else {
2419:     if (set) *set = PETSC_FALSE;
2420:   }
2421:   PetscFunctionReturn(PETSC_SUCCESS);
2422: }

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

2428:   Not Collective

2430:   Input Parameters:
2431: + options - options database, use `NULL` for default global database
2432: . pre     - the string to prepend to the name or `NULL`
2433: - name    - the option one is seeking

2435:   Output Parameters:
2436: + ivalue - the ternary logical value to return
2437: - set    - `PETSC_TRUE`  if found, else `PETSC_FALSE`

2439:   Level: beginner

2441:   Notes:
2442:   TRUE, true, YES, yes, ON, on, nostring and 1 all translate to `PETSC_BOOL3_TRUE`
2443:   FALSE, false, NO, no, OFF, off and 0 all translate to `PETSC_BOOL3_FALSE`
2444:   UNKNOWN, unknown, AUTO and auto all translate to `PETSC_BOOL3_UNKNOWN`

2446:   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`
2447:   is equivalent to `-requested_bool3 true`

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

2452: .seealso: `PetscOptionsGetBool()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2453:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetInt()`, `PetscOptionsBool()`,
2454:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2455:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2456:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2457:           `PetscOptionsFList()`, `PetscOptionsEList()`
2458: @*/
2459: PetscErrorCode PetscOptionsGetBool3(PetscOptions options, const char pre[], const char name[], PetscBool3 *ivalue, PetscBool *set)
2460: {
2461:   const char *value;
2462:   PetscBool   flag;

2464:   PetscFunctionBegin;
2465:   PetscAssertPointer(name, 3);
2466:   if (ivalue) PetscAssertPointer(ivalue, 4);
2467:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2468:   if (flag) { // found the option
2469:     PetscBool isAUTO = PETSC_FALSE, isUNKNOWN = PETSC_FALSE;

2471:     if (set) *set = PETSC_TRUE;
2472:     PetscCall(PetscStrcasecmp("AUTO", value, &isAUTO));                    // auto or AUTO
2473:     if (!isAUTO) PetscCall(PetscStrcasecmp("UNKNOWN", value, &isUNKNOWN)); // unknown or UNKNOWN
2474:     if (isAUTO || isUNKNOWN) {
2475:       if (ivalue) *ivalue = PETSC_BOOL3_UNKNOWN;
2476:     } else { // handle boolean values (if no value is given, it returns true)
2477:       PetscCall(PetscOptionsStringToBool(value, &flag));
2478:       if (ivalue) *ivalue = PetscBoolToBool3(flag);
2479:     }
2480:   } else {
2481:     if (set) *set = PETSC_FALSE;
2482:   }
2483:   PetscFunctionReturn(PETSC_SUCCESS);
2484: }

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

2489:   Not Collective

2491:   Input Parameters:
2492: + options - options database, use `NULL` for default global database
2493: . pre     - the string to prepend to the name or `NULL`
2494: . opt     - option name
2495: . list    - the possible choices (one of these must be selected, anything else is invalid)
2496: - ntext   - number of choices

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

2502:   Level: intermediate

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

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

2510: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2511:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2512:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2513:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2514:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2515:           `PetscOptionsFList()`, `PetscOptionsEList()`
2516: @*/
2517: PetscErrorCode PetscOptionsGetEList(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscInt ntext, PetscInt *value, PetscBool *set)
2518: {
2519:   size_t    alen, len = 0, tlen = 0;
2520:   char     *svalue;
2521:   PetscBool aset, flg = PETSC_FALSE;
2522:   PetscInt  i;

2524:   PetscFunctionBegin;
2525:   PetscAssertPointer(opt, 3);
2526:   for (i = 0; i < ntext; i++) {
2527:     PetscCall(PetscStrlen(list[i], &alen));
2528:     if (alen > len) len = alen;
2529:     tlen += len + 1;
2530:   }
2531:   len += 5; /* a little extra space for user mistypes */
2532:   PetscCall(PetscMalloc1(len, &svalue));
2533:   PetscCall(PetscOptionsGetString(options, pre, opt, svalue, len, &aset));
2534:   if (aset) {
2535:     PetscCall(PetscEListFind(ntext, list, svalue, value, &flg));
2536:     if (!flg) {
2537:       char *avail;

2539:       PetscCall(PetscMalloc1(tlen, &avail));
2540:       avail[0] = '\0';
2541:       for (i = 0; i < ntext; i++) {
2542:         PetscCall(PetscStrlcat(avail, list[i], tlen));
2543:         PetscCall(PetscStrlcat(avail, " ", tlen));
2544:       }
2545:       PetscCall(PetscStrtolower(avail));
2546:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown option %s for -%s%s. Available options: %s", svalue, pre ? pre : "", opt + 1, avail);
2547:     }
2548:     if (set) *set = PETSC_TRUE;
2549:   } else if (set) *set = PETSC_FALSE;
2550:   PetscCall(PetscFree(svalue));
2551:   PetscFunctionReturn(PETSC_SUCCESS);
2552: }

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

2557:   Not Collective

2559:   Input Parameters:
2560: + options - options database, use `NULL` for default global database
2561: . pre     - option prefix or `NULL`
2562: . opt     - option name
2563: - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2565:   Output Parameters:
2566: + value - the value to return
2567: - set   - `PETSC_TRUE` if found, else `PETSC_FALSE`

2569:   Level: beginner

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

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

2577: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2578:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2579:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2580:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2581:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2582:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2583:           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2584: @*/
2585: PetscErrorCode PetscOptionsGetEnum(PetscOptions options, const char pre[], const char opt[], const char *const *list, PetscEnum *value, PetscBool *set)
2586: {
2587:   PetscInt  ntext = 0, tval;
2588:   PetscBool fset;

2590:   PetscFunctionBegin;
2591:   PetscAssertPointer(opt, 3);
2592:   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");
2593:   PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix");
2594:   ntext -= 3;
2595:   PetscCall(PetscOptionsGetEList(options, pre, opt, list, ntext, &tval, &fset));
2596:   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2597:   if (fset) *value = (PetscEnum)tval;
2598:   if (set) *set = fset;
2599:   PetscFunctionReturn(PETSC_SUCCESS);
2600: }

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

2605:   Not Collective

2607:   Input Parameters:
2608: + options - options database, use `NULL` for default global database
2609: . pre     - the string to prepend to the name or `NULL`
2610: - name    - the option one is seeking

2612:   Output Parameters:
2613: + ivalue - the integer value to return
2614: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2616:   Level: beginner

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

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

2624:   Accepts the deprecated value `default`.

2626: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2627:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2628:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2629:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2630:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2631:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2632:           `PetscOptionsFList()`, `PetscOptionsEList()`
2633: @*/
2634: PetscErrorCode PetscOptionsGetInt(PetscOptions options, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
2635: {
2636:   const char *value;
2637:   PetscBool   flag;

2639:   PetscFunctionBegin;
2640:   PetscAssertPointer(name, 3);
2641:   PetscAssertPointer(ivalue, 4);
2642:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2643:   if (flag) {
2644:     if (!value) {
2645:       if (set) *set = PETSC_FALSE;
2646:     } else {
2647:       if (set) *set = PETSC_TRUE;
2648:       PetscCall(PetscOptionsStringToInt(value, ivalue));
2649:     }
2650:   } else {
2651:     if (set) *set = PETSC_FALSE;
2652:   }
2653:   PetscFunctionReturn(PETSC_SUCCESS);
2654: }

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

2659:   Not Collective

2661:   Input Parameters:
2662: + options - options database, use `NULL` for default global database
2663: . pre     - the string to prepend to the name or `NULL`
2664: - name    - the option one is seeking

2666:   Output Parameters:
2667: + ivalue - the MPI integer value to return
2668: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2670:   Level: beginner

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

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

2678:   Accepts the deprecated value `default`.

2680: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
2681:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2682:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
2683:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2684:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2685:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2686:           `PetscOptionsFList()`, `PetscOptionsEList()`
2687: @*/
2688: PetscErrorCode PetscOptionsGetMPIInt(PetscOptions options, const char pre[], const char name[], PetscMPIInt *ivalue, PetscBool *set)
2689: {
2690:   PetscInt  value;
2691:   PetscBool flag;

2693:   PetscFunctionBegin;
2694:   PetscCall(PetscOptionsGetInt(options, pre, name, &value, &flag));
2695:   if (flag) PetscCall(PetscMPIIntCast(value, ivalue));
2696:   if (set) *set = flag;
2697:   PetscFunctionReturn(PETSC_SUCCESS);
2698: }

2700: /*@C
2701:   PetscOptionsGetReal - Gets the double precision value for a particular
2702:   option in the database.

2704:   Not Collective

2706:   Input Parameters:
2707: + options - options database, use `NULL` for default global database
2708: . pre     - string to prepend to each name or `NULL`
2709: - name    - the option one is seeking

2711:   Output Parameters:
2712: + dvalue - the double value to return
2713: - set    - `PETSC_TRUE` if found, `PETSC_FALSE` if not found

2715:   Level: beginner

2717:   Notes:
2718:   Accepts the special values `determine`, `decide` and `unlimited`.

2720:   Accepts the deprecated value `default`

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

2725: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2726:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2727:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2728:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2729:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2730:           `PetscOptionsFList()`, `PetscOptionsEList()`
2731: @*/
2732: PetscErrorCode PetscOptionsGetReal(PetscOptions options, const char pre[], const char name[], PetscReal *dvalue, PetscBool *set)
2733: {
2734:   const char *value;
2735:   PetscBool   flag;

2737:   PetscFunctionBegin;
2738:   PetscAssertPointer(name, 3);
2739:   PetscAssertPointer(dvalue, 4);
2740:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2741:   if (flag) {
2742:     if (!value) {
2743:       if (set) *set = PETSC_FALSE;
2744:     } else {
2745:       if (set) *set = PETSC_TRUE;
2746:       PetscCall(PetscOptionsStringToReal(value, dvalue));
2747:     }
2748:   } else {
2749:     if (set) *set = PETSC_FALSE;
2750:   }
2751:   PetscFunctionReturn(PETSC_SUCCESS);
2752: }

2754: /*@C
2755:   PetscOptionsGetScalar - Gets the scalar value for a particular
2756:   option in the database.

2758:   Not Collective

2760:   Input Parameters:
2761: + options - options database, use `NULL` for default global database
2762: . pre     - string to prepend to each name or `NULL`
2763: - name    - the option one is seeking

2765:   Output Parameters:
2766: + dvalue - the scalar value to return
2767: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2769:   Level: beginner

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

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

2778: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2779:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2780:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2781:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2782:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2783:           `PetscOptionsFList()`, `PetscOptionsEList()`
2784: @*/
2785: PetscErrorCode PetscOptionsGetScalar(PetscOptions options, const char pre[], const char name[], PetscScalar *dvalue, PetscBool *set)
2786: {
2787:   const char *value;
2788:   PetscBool   flag;

2790:   PetscFunctionBegin;
2791:   PetscAssertPointer(name, 3);
2792:   PetscAssertPointer(dvalue, 4);
2793:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2794:   if (flag) {
2795:     if (!value) {
2796:       if (set) *set = PETSC_FALSE;
2797:     } else {
2798: #if !defined(PETSC_USE_COMPLEX)
2799:       PetscCall(PetscOptionsStringToReal(value, dvalue));
2800: #else
2801:       PetscCall(PetscOptionsStringToScalar(value, dvalue));
2802: #endif
2803:       if (set) *set = PETSC_TRUE;
2804:     }
2805:   } else { /* flag */
2806:     if (set) *set = PETSC_FALSE;
2807:   }
2808:   PetscFunctionReturn(PETSC_SUCCESS);
2809: }

2811: /*@C
2812:   PetscOptionsGetString - Gets the string value for a particular option in
2813:   the database.

2815:   Not Collective

2817:   Input Parameters:
2818: + options - options database, use `NULL` for default global database
2819: . pre     - string to prepend to name or `NULL`
2820: . name    - the option one is seeking
2821: - len     - maximum length of the string including null termination

2823:   Output Parameters:
2824: + string - location to copy string
2825: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2827:   Level: beginner

2829:   Note:
2830:   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`

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

2835:   Fortran Notes:
2836:   The Fortran interface is slightly different from the C/C++
2837:   interface (len is not used).  Sample usage in Fortran follows
2838: .vb
2839:       character *20    string
2840:       PetscErrorCode   ierr
2841:       PetscBool        set
2842:       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2843: .ve

2845: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
2846:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2847:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2848:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2849:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2850:           `PetscOptionsFList()`, `PetscOptionsEList()`
2851: @*/
2852: PetscErrorCode PetscOptionsGetString(PetscOptions options, const char pre[], const char name[], char string[], size_t len, PetscBool *set)
2853: {
2854:   const char *value;
2855:   PetscBool   flag;

2857:   PetscFunctionBegin;
2858:   PetscAssertPointer(name, 3);
2859:   PetscAssertPointer(string, 4);
2860:   PetscCall(PetscOptionsFindPair(options, pre, name, &value, &flag));
2861:   if (!flag) {
2862:     if (set) *set = PETSC_FALSE;
2863:   } else {
2864:     if (set) *set = PETSC_TRUE;
2865:     if (value) PetscCall(PetscStrncpy(string, value, len));
2866:     else PetscCall(PetscArrayzero(string, len));
2867:   }
2868:   PetscFunctionReturn(PETSC_SUCCESS);
2869: }

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

2875:   Not Collective

2877:   Input Parameters:
2878: + options - options database, use `NULL` for default global database
2879: . pre     - string to prepend to each name or `NULL`
2880: - name    - the option one is seeking

2882:   Output Parameters:
2883: + dvalue - the Boolean values to return
2884: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2885: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2887:   Level: beginner

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

2892: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
2893:           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
2894:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
2895:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
2896:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2897:           `PetscOptionsFList()`, `PetscOptionsEList()`
2898: @*/
2899: PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options, const char pre[], const char name[], PetscBool dvalue[], PetscInt *nmax, PetscBool *set)
2900: {
2901:   const char *svalue;
2902:   char       *value;
2903:   PetscInt    n = 0;
2904:   PetscBool   flag;
2905:   PetscToken  token;

2907:   PetscFunctionBegin;
2908:   PetscAssertPointer(name, 3);
2909:   PetscAssertPointer(dvalue, 4);
2910:   PetscAssertPointer(nmax, 5);

2912:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2913:   if (!flag || !svalue) {
2914:     if (set) *set = PETSC_FALSE;
2915:     *nmax = 0;
2916:     PetscFunctionReturn(PETSC_SUCCESS);
2917:   }
2918:   if (set) *set = PETSC_TRUE;
2919:   PetscCall(PetscTokenCreate(svalue, ',', &token));
2920:   PetscCall(PetscTokenFind(token, &value));
2921:   while (value && n < *nmax) {
2922:     PetscCall(PetscOptionsStringToBool(value, dvalue));
2923:     PetscCall(PetscTokenFind(token, &value));
2924:     dvalue++;
2925:     n++;
2926:   }
2927:   PetscCall(PetscTokenDestroy(&token));
2928:   *nmax = n;
2929:   PetscFunctionReturn(PETSC_SUCCESS);
2930: }

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

2935:   Not Collective

2937:   Input Parameters:
2938: + options - options database, use `NULL` for default global database
2939: . pre     - option prefix or `NULL`
2940: . name    - option name
2941: - list    - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2943:   Output Parameters:
2944: + ivalue - the  enum values to return
2945: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
2946: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

2948:   Level: beginner

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

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

2955: .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`,
2956:           `PetscOptionsGetEnum()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
2957:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsName()`,
2958:           `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, `PetscOptionsStringArray()`, `PetscOptionsRealArray()`,
2959:           `PetscOptionsScalar()`, `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
2960:           `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsGetEList()`, `PetscOptionsEnum()`
2961: @*/
2962: PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options, const char pre[], const char name[], const char *const *list, PetscEnum ivalue[], PetscInt *nmax, PetscBool *set)
2963: {
2964:   const char *svalue;
2965:   char       *value;
2966:   PetscInt    n = 0;
2967:   PetscEnum   evalue;
2968:   PetscBool   flag;
2969:   PetscToken  token;

2971:   PetscFunctionBegin;
2972:   PetscAssertPointer(name, 3);
2973:   PetscAssertPointer(list, 4);
2974:   PetscAssertPointer(ivalue, 5);
2975:   PetscAssertPointer(nmax, 6);

2977:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
2978:   if (!flag || !svalue) {
2979:     if (set) *set = PETSC_FALSE;
2980:     *nmax = 0;
2981:     PetscFunctionReturn(PETSC_SUCCESS);
2982:   }
2983:   if (set) *set = PETSC_TRUE;
2984:   PetscCall(PetscTokenCreate(svalue, ',', &token));
2985:   PetscCall(PetscTokenFind(token, &value));
2986:   while (value && n < *nmax) {
2987:     PetscCall(PetscEnumFind(list, value, &evalue, &flag));
2988:     PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown enum value '%s' for -%s%s", svalue, pre ? pre : "", name + 1);
2989:     ivalue[n++] = evalue;
2990:     PetscCall(PetscTokenFind(token, &value));
2991:   }
2992:   PetscCall(PetscTokenDestroy(&token));
2993:   *nmax = n;
2994:   PetscFunctionReturn(PETSC_SUCCESS);
2995: }

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

3000:   Not Collective

3002:   Input Parameters:
3003: + options - options database, use `NULL` for default global database
3004: . pre     - string to prepend to each name or `NULL`
3005: - name    - the option one is seeking

3007:   Output Parameters:
3008: + ivalue - the integer values to return
3009: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3010: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3012:   Level: beginner

3014:   Notes:
3015:   The array can be passed as
3016: +  a comma separated list -                                 0,1,2,3,4,5,6,7
3017: .  a range (start\-end+1) -                                 0-8
3018: .  a range with given increment (start\-end+1:inc) -        0-7:2
3019: -  a combination of values and ranges separated by commas - 0,1-8,8-15:2

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

3023: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3024:           `PetscOptionsGetString()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3025:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3026:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3027:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3028:           `PetscOptionsFList()`, `PetscOptionsEList()`
3029: @*/
3030: PetscErrorCode PetscOptionsGetIntArray(PetscOptions options, const char pre[], const char name[], PetscInt ivalue[], PetscInt *nmax, PetscBool *set)
3031: {
3032:   const char *svalue;
3033:   char       *value;
3034:   PetscInt    n = 0, i, j, start, end, inc, nvalues;
3035:   size_t      len;
3036:   PetscBool   flag, foundrange;
3037:   PetscToken  token;

3039:   PetscFunctionBegin;
3040:   PetscAssertPointer(name, 3);
3041:   PetscAssertPointer(ivalue, 4);
3042:   PetscAssertPointer(nmax, 5);

3044:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3045:   if (!flag || !svalue) {
3046:     if (set) *set = PETSC_FALSE;
3047:     *nmax = 0;
3048:     PetscFunctionReturn(PETSC_SUCCESS);
3049:   }
3050:   if (set) *set = PETSC_TRUE;
3051:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3052:   PetscCall(PetscTokenFind(token, &value));
3053:   while (value && n < *nmax) {
3054:     /* look for form  d-D where d and D are integers */
3055:     foundrange = PETSC_FALSE;
3056:     PetscCall(PetscStrlen(value, &len));
3057:     if (value[0] == '-') i = 2;
3058:     else i = 1;
3059:     for (; i < (int)len; i++) {
3060:       if (value[i] == '-') {
3061:         PetscCheck(i != (int)len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value);
3062:         value[i] = 0;

3064:         PetscCall(PetscOptionsStringToInt(value, &start));
3065:         inc = 1;
3066:         j   = i + 1;
3067:         for (; j < (int)len; j++) {
3068:           if (value[j] == ':') {
3069:             value[j] = 0;

3071:             PetscCall(PetscOptionsStringToInt(value + j + 1, &inc));
3072:             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);
3073:             break;
3074:           }
3075:         }
3076:         PetscCall(PetscOptionsStringToInt(value + i + 1, &end));
3077:         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);
3078:         nvalues = (end - start) / inc + (end - start) % inc;
3079:         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);
3080:         for (; start < end; start += inc) {
3081:           *ivalue = start;
3082:           ivalue++;
3083:           n++;
3084:         }
3085:         foundrange = PETSC_TRUE;
3086:         break;
3087:       }
3088:     }
3089:     if (!foundrange) {
3090:       PetscCall(PetscOptionsStringToInt(value, ivalue));
3091:       ivalue++;
3092:       n++;
3093:     }
3094:     PetscCall(PetscTokenFind(token, &value));
3095:   }
3096:   PetscCall(PetscTokenDestroy(&token));
3097:   *nmax = n;
3098:   PetscFunctionReturn(PETSC_SUCCESS);
3099: }

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

3105:   Not Collective

3107:   Input Parameters:
3108: + options - options database, use `NULL` for default global database
3109: . pre     - string to prepend to each name or `NULL`
3110: - name    - the option one is seeking

3112:   Output Parameters:
3113: + dvalue - the double values to return
3114: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3115: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3117:   Level: beginner

3119: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3120:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3121:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3122:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3123:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3124:           `PetscOptionsFList()`, `PetscOptionsEList()`
3125: @*/
3126: PetscErrorCode PetscOptionsGetRealArray(PetscOptions options, const char pre[], const char name[], PetscReal dvalue[], PetscInt *nmax, PetscBool *set)
3127: {
3128:   const char *svalue;
3129:   char       *value;
3130:   PetscInt    n = 0;
3131:   PetscBool   flag;
3132:   PetscToken  token;

3134:   PetscFunctionBegin;
3135:   PetscAssertPointer(name, 3);
3136:   PetscAssertPointer(dvalue, 4);
3137:   PetscAssertPointer(nmax, 5);

3139:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3140:   if (!flag || !svalue) {
3141:     if (set) *set = PETSC_FALSE;
3142:     *nmax = 0;
3143:     PetscFunctionReturn(PETSC_SUCCESS);
3144:   }
3145:   if (set) *set = PETSC_TRUE;
3146:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3147:   PetscCall(PetscTokenFind(token, &value));
3148:   while (value && n < *nmax) {
3149:     PetscCall(PetscOptionsStringToReal(value, dvalue++));
3150:     PetscCall(PetscTokenFind(token, &value));
3151:     n++;
3152:   }
3153:   PetscCall(PetscTokenDestroy(&token));
3154:   *nmax = n;
3155:   PetscFunctionReturn(PETSC_SUCCESS);
3156: }

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

3162:   Not Collective

3164:   Input Parameters:
3165: + options - options database, use `NULL` for default global database
3166: . pre     - string to prepend to each name or `NULL`
3167: - name    - the option one is seeking

3169:   Output Parameters:
3170: + dvalue - the scalar values to return
3171: . nmax   - On input maximum number of values to retrieve, on output the actual number of values retrieved
3172: - set    - `PETSC_TRUE` if found, else `PETSC_FALSE`

3174:   Level: beginner

3176: .seealso: `PetscOptionsGetInt()`, `PetscOptionsHasName()`,
3177:           `PetscOptionsGetString()`, `PetscOptionsGetIntArray()`, `PetscOptionsBool()`,
3178:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3179:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3180:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3181:           `PetscOptionsFList()`, `PetscOptionsEList()`
3182: @*/
3183: PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options, const char pre[], const char name[], PetscScalar dvalue[], PetscInt *nmax, PetscBool *set)
3184: {
3185:   const char *svalue;
3186:   char       *value;
3187:   PetscInt    n = 0;
3188:   PetscBool   flag;
3189:   PetscToken  token;

3191:   PetscFunctionBegin;
3192:   PetscAssertPointer(name, 3);
3193:   PetscAssertPointer(dvalue, 4);
3194:   PetscAssertPointer(nmax, 5);

3196:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3197:   if (!flag || !svalue) {
3198:     if (set) *set = PETSC_FALSE;
3199:     *nmax = 0;
3200:     PetscFunctionReturn(PETSC_SUCCESS);
3201:   }
3202:   if (set) *set = PETSC_TRUE;
3203:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3204:   PetscCall(PetscTokenFind(token, &value));
3205:   while (value && n < *nmax) {
3206:     PetscCall(PetscOptionsStringToScalar(value, dvalue++));
3207:     PetscCall(PetscTokenFind(token, &value));
3208:     n++;
3209:   }
3210:   PetscCall(PetscTokenDestroy(&token));
3211:   *nmax = n;
3212:   PetscFunctionReturn(PETSC_SUCCESS);
3213: }

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

3219:   Not Collective; No Fortran Support

3221:   Input Parameters:
3222: + options - options database, use `NULL` for default global database
3223: . pre     - string to prepend to name or `NULL`
3224: - name    - the option one is seeking

3226:   Output Parameters:
3227: + strings - location to copy strings
3228: . nmax    - On input maximum number of strings, on output the actual number of strings found
3229: - set     - `PETSC_TRUE` if found, else `PETSC_FALSE`

3231:   Level: beginner

3233:   Notes:
3234:   The `nmax` parameter is used for both input and output.

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

3239:   The user is responsible for deallocating the strings that are
3240:   returned.

3242: .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`,
3243:           `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
3244:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
3245:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
3246:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
3247:           `PetscOptionsFList()`, `PetscOptionsEList()`
3248: @*/
3249: PetscErrorCode PetscOptionsGetStringArray(PetscOptions options, const char pre[], const char name[], char *strings[], PetscInt *nmax, PetscBool *set)
3250: {
3251:   const char *svalue;
3252:   char       *value;
3253:   PetscInt    n = 0;
3254:   PetscBool   flag;
3255:   PetscToken  token;

3257:   PetscFunctionBegin;
3258:   PetscAssertPointer(name, 3);
3259:   PetscAssertPointer(strings, 4);
3260:   PetscAssertPointer(nmax, 5);

3262:   PetscCall(PetscOptionsFindPair(options, pre, name, &svalue, &flag));
3263:   if (!flag || !svalue) {
3264:     if (set) *set = PETSC_FALSE;
3265:     *nmax = 0;
3266:     PetscFunctionReturn(PETSC_SUCCESS);
3267:   }
3268:   if (set) *set = PETSC_TRUE;
3269:   PetscCall(PetscTokenCreate(svalue, ',', &token));
3270:   PetscCall(PetscTokenFind(token, &value));
3271:   while (value && n < *nmax) {
3272:     PetscCall(PetscStrallocpy(value, &strings[n]));
3273:     PetscCall(PetscTokenFind(token, &value));
3274:     n++;
3275:   }
3276:   PetscCall(PetscTokenDestroy(&token));
3277:   *nmax = n;
3278:   PetscFunctionReturn(PETSC_SUCCESS);
3279: }

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

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

3286:   Logically Collective

3288:   Input Parameters:
3289: + PetscOptionsObject - string to prepend to name or `NULL`
3290: . oldname            - the old, deprecated option
3291: . newname            - the new option, or `NULL` if option is purely removed
3292: . version            - a string describing the version of first deprecation, e.g. "3.9"
3293: - info               - additional information string, or `NULL`.

3295:   Options Database Key:
3296: . -options_suppress_deprecated_warnings - do not print deprecation warnings

3298:   Level: developer

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

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

3307:   Must be called between `PetscOptionsBegin()` (or `PetscObjectOptionsBegin()`) and `PetscOptionsEnd()`.
3308:   Only the process of rank zero that owns the `PetscOptionsItems` are argument (managed by `PetscOptionsBegin()` or
3309:   `PetscObjectOptionsBegin()` prints the information
3310:   If newname is provided, the old option is replaced. Otherwise, it remains
3311:   in the options database.
3312:   If an option is not replaced, the info argument should be used to advise the user
3313:   on how to proceed.
3314:   There is a limit on the length of the warning printed, so very long strings
3315:   provided as info may be truncated.

3317: .seealso: `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsScalar()`, `PetscOptionsBool()`, `PetscOptionsString()`, `PetscOptionsSetValue()`
3318: @*/
3319: PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject, const char oldname[], const char newname[], const char version[], const char info[])
3320: {
3321:   PetscBool         found, quiet;
3322:   const char       *value;
3323:   const char *const quietopt = "-options_suppress_deprecated_warnings";
3324:   char              msg[4096];
3325:   char             *prefix  = NULL;
3326:   PetscOptions      options = NULL;
3327:   MPI_Comm          comm    = PETSC_COMM_SELF;

3329:   PetscFunctionBegin;
3330:   PetscAssertPointer(oldname, 2);
3331:   PetscAssertPointer(version, 4);
3332:   if (PetscOptionsObject) {
3333:     prefix  = PetscOptionsObject->prefix;
3334:     options = PetscOptionsObject->options;
3335:     comm    = PetscOptionsObject->comm;
3336:   }
3337:   PetscCall(PetscOptionsFindPair(options, prefix, oldname, &value, &found));
3338:   if (found) {
3339:     if (newname) {
3340:       if (prefix) PetscCall(PetscOptionsPrefixPush(options, prefix));
3341:       PetscCall(PetscOptionsSetValue(options, newname, value));
3342:       if (prefix) PetscCall(PetscOptionsPrefixPop(options));
3343:       PetscCall(PetscOptionsClearValue(options, oldname));
3344:     }
3345:     quiet = PETSC_FALSE;
3346:     PetscCall(PetscOptionsGetBool(options, NULL, quietopt, &quiet, NULL));
3347:     if (!quiet) {
3348:       PetscCall(PetscStrncpy(msg, "** PETSc DEPRECATION WARNING ** : the option -", sizeof(msg)));
3349:       PetscCall(PetscStrlcat(msg, prefix, sizeof(msg)));
3350:       PetscCall(PetscStrlcat(msg, oldname + 1, sizeof(msg)));
3351:       PetscCall(PetscStrlcat(msg, " is deprecated as of version ", sizeof(msg)));
3352:       PetscCall(PetscStrlcat(msg, version, sizeof(msg)));
3353:       PetscCall(PetscStrlcat(msg, " and will be removed in a future release.\n", sizeof(msg)));
3354:       if (newname) {
3355:         PetscCall(PetscStrlcat(msg, "   Use the option -", sizeof(msg)));
3356:         PetscCall(PetscStrlcat(msg, prefix, sizeof(msg)));
3357:         PetscCall(PetscStrlcat(msg, newname + 1, sizeof(msg)));
3358:         PetscCall(PetscStrlcat(msg, " instead.", sizeof(msg)));
3359:       }
3360:       if (info) {
3361:         PetscCall(PetscStrlcat(msg, " ", sizeof(msg)));
3362:         PetscCall(PetscStrlcat(msg, info, sizeof(msg)));
3363:       }
3364:       PetscCall(PetscStrlcat(msg, " (Silence this warning with ", sizeof(msg)));
3365:       PetscCall(PetscStrlcat(msg, quietopt, sizeof(msg)));
3366:       PetscCall(PetscStrlcat(msg, ")\n", sizeof(msg)));
3367:       PetscCall(PetscPrintf(comm, "%s", msg));
3368:     }
3369:   }
3370:   PetscFunctionReturn(PETSC_SUCCESS);
3371: }