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