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