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