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