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