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