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