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