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