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