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