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