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