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