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