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