xref: /petsc/src/sys/objects/options.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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 /*@
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 .  viewer - must be an PETSCVIEWERASCII viewer
757 
758    Options Database Key:
759 .  -options_table - Activates PetscOptionsView() within PetscFinalize()
760 
761    Level: advanced
762 
763    Concepts: options database^printing
764 
765 .seealso: PetscOptionsAllUsed()
766 @*/
767 PetscErrorCode  PetscOptionsView(PetscOptions options,PetscViewer viewer)
768 {
769   PetscErrorCode ierr;
770   PetscInt       i;
771   PetscBool      isascii;
772 
773   PetscFunctionBegin;
774   options = options ? options : defaultoptions;
775   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
776   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
777   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Only supports ASCII viewer");
778 
779   if (options->N) {
780     ierr = PetscViewerASCIIPrintf(viewer,"#PETSc Option Table entries:\n");CHKERRQ(ierr);
781   } else {
782     ierr = PetscViewerASCIIPrintf(viewer,"#No PETSc Option Table entries\n");CHKERRQ(ierr);
783   }
784   for (i=0; i<options->N; i++) {
785     if (options->values[i]) {
786       ierr = PetscViewerASCIIPrintf(viewer,"-%s %s\n",options->names[i],options->values[i]);CHKERRQ(ierr);
787     } else {
788       ierr = PetscViewerASCIIPrintf(viewer,"-%s\n",options->names[i]);CHKERRQ(ierr);
789     }
790   }
791   if (options->N) {
792     ierr = PetscViewerASCIIPrintf(viewer,"#End of PETSc Option Table entries\n");CHKERRQ(ierr);
793   }
794   PetscFunctionReturn(0);
795 }
796 
797 /*
798    Called by error handlers to print options used in run
799 */
800 PetscErrorCode  PetscOptionsViewError(void)
801 {
802   PetscInt       i;
803   PetscOptions   options = defaultoptions;
804 
805   PetscFunctionBegin;
806   if (options->N) {
807     (*PetscErrorPrintf)("PETSc Option Table entries:\n");
808   } else {
809     (*PetscErrorPrintf)("No PETSc Option Table entries\n");
810   }
811   for (i=0; i<options->N; i++) {
812     if (options->values[i]) {
813       (*PetscErrorPrintf)("-%s %s\n",options->names[i],options->values[i]);
814     } else {
815       (*PetscErrorPrintf)("-%s\n",options->names[i]);
816     }
817   }
818   PetscFunctionReturn(0);
819 }
820 
821 /*@C
822    PetscOptionsGetAll - Lists all the options the program was run with in a single string.
823 
824    Not Collective
825 
826    Input Paramter:
827 .  options - the options database, use NULL for the default global database
828 
829    Output Parameter:
830 .  copts - pointer where string pointer is stored
831 
832    Notes: the array and each entry in the array should be freed with PetscFree()
833 
834    Level: advanced
835 
836    Concepts: options database^listing
837 
838 .seealso: PetscOptionsAllUsed(), PetscOptionsView()
839 @*/
840 PetscErrorCode  PetscOptionsGetAll(PetscOptions options,char *copts[])
841 {
842   PetscErrorCode ierr;
843   PetscInt       i;
844   size_t         len       = 1,lent = 0;
845   char           *coptions = NULL;
846 
847   PetscFunctionBegin;
848   options = options ? options : defaultoptions;
849 
850   /* count the length of the required string */
851   for (i=0; i<options->N; i++) {
852     ierr = PetscStrlen(options->names[i],&lent);CHKERRQ(ierr);
853     len += 2 + lent;
854     if (options->values[i]) {
855       ierr = PetscStrlen(options->values[i],&lent);CHKERRQ(ierr);
856       len += 1 + lent;
857     }
858   }
859   ierr = PetscMalloc1(len,&coptions);CHKERRQ(ierr);
860   coptions[0] = 0;
861   for (i=0; i<options->N; i++) {
862     ierr = PetscStrcat(coptions,"-");CHKERRQ(ierr);
863     ierr = PetscStrcat(coptions,options->names[i]);CHKERRQ(ierr);
864     ierr = PetscStrcat(coptions," ");CHKERRQ(ierr);
865     if (options->values[i]) {
866       ierr = PetscStrcat(coptions,options->values[i]);CHKERRQ(ierr);
867       ierr = PetscStrcat(coptions," ");CHKERRQ(ierr);
868     }
869   }
870   *copts = coptions;
871   PetscFunctionReturn(0);
872 }
873 
874 /*@C
875    PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.
876 
877    Not Collective, but prefix will only be applied on calling ranks
878 
879    Input Parameter:
880 +  options - options database, or NULL for the default global database
881 -  prefix - The string to append to the existing prefix
882 
883    Options Database Keys:
884  +   -prefix_push <some_prefix_> - push the given prefix
885  -   -prefix_pop - pop the last prefix
886 
887    Notes:
888    It is common to use this in conjunction with -options_file as in
889 
890  $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
891 
892    where the files no longer require all options to be prefixed with -system2_.
893 
894 Level: advanced
895 
896 .seealso: PetscOptionsPrefixPop()
897 @*/
898 PetscErrorCode  PetscOptionsPrefixPush(PetscOptions options,const char prefix[])
899 {
900   PetscErrorCode ierr;
901   size_t         n;
902   PetscInt       start;
903   char           buf[2048];
904   PetscBool      key;
905 
906   PetscFunctionBegin;
907   PetscValidCharPointer(prefix,1);
908   options = options ? options : defaultoptions;
909   /* Want to check validity of the key using PetscOptionsValidKey(), which requires that the first character is a '-' */
910   buf[0] = '-';
911   ierr = PetscStrncpy(buf+1,prefix,sizeof(buf) - 1);CHKERRQ(ierr);
912   buf[sizeof(buf) - 1] = 0;
913   ierr = PetscOptionsValidKey(buf,&key);CHKERRQ(ierr);
914   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);
915 
916   if (!options) {
917     ierr = PetscOptionsInsert(NULL,0,0,0);CHKERRQ(ierr);
918     options = defaultoptions;
919   }
920   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);
921   start = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
922   ierr = PetscStrlen(prefix,&n);CHKERRQ(ierr);
923   if (n+1 > sizeof(options->prefix)-start) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Maximum prefix length %d exceeded",sizeof(options->prefix));
924   ierr = PetscMemcpy(options->prefix+start,prefix,n+1);CHKERRQ(ierr);
925   options->prefixstack[options->prefixind++] = start+n;
926   PetscFunctionReturn(0);
927 }
928 
929 /*@C
930    PetscOptionsPrefixPop - Remove the latest options prefix, see PetscOptionsPrefixPush() for details
931 
932    Not  Collective, but prefix will only be popped on calling ranks
933 
934   Input Parameters:
935 .  options - options database, or NULL for the default global database
936 
937    Level: advanced
938 
939 .seealso: PetscOptionsPrefixPush()
940 @*/
941 PetscErrorCode  PetscOptionsPrefixPop(PetscOptions options)
942 {
943   PetscInt offset;
944 
945   PetscFunctionBegin;
946   options = options ? options : defaultoptions;
947   if (options->prefixind < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"More prefixes popped than pushed");
948   options->prefixind--;
949   offset = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
950   options->prefix[offset] = 0;
951   PetscFunctionReturn(0);
952 }
953 
954 /*@C
955     PetscOptionsClear - Removes all options form the database leaving it empty.
956 
957   Input Parameters:
958 .  options - options database, use NULL for the default global database
959 
960    Level: developer
961 
962 .seealso: PetscOptionsInsert()
963 @*/
964 PetscErrorCode  PetscOptionsClear(PetscOptions options)
965 {
966   PetscInt     i;
967 
968   PetscFunctionBegin;
969   options = options ? options : defaultoptions;
970   for (i=0; i<options->N; i++) {
971     if (options->names[i])  free(options->names[i]);
972     if (options->values[i]) free(options->values[i]);
973   }
974   for (i=0; i<options->Naliases; i++) {
975     free(options->aliases1[i]);
976     free(options->aliases2[i]);
977   }
978   options->prefix[0] = 0;
979   options->prefixind = 0;
980   options->N         = 0;
981   options->Naliases  = 0;
982   PetscFunctionReturn(0);
983 }
984 
985 /*@
986     PetscObjectSetPrintedOptions - indicate to an object that it should behave as if it has already printed the help for its options
987 
988   Input Parameters:
989 .  obj  - the PetscObject
990 
991    Level: developer
992 
993    Developer Notes: This is used, for example to prevent sequential objects that are created from a parallel object; such as the KSP created by
994     PCBJACOBI from all printing the same help messages to the screen
995 
996 .seealso: PetscOptionsInsert()
997 @*/
998 PetscErrorCode  PetscObjectSetPrintedOptions(PetscObject obj)
999 {
1000   PetscFunctionBegin;
1001   obj->optionsprinted = PETSC_TRUE;
1002   PetscFunctionReturn(0);
1003 }
1004 
1005 /*@
1006     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.
1007 
1008   Input Parameters:
1009 +  pobj - the parent object
1010 -  obj  - the PetscObject
1011 
1012    Level: developer
1013 
1014    Developer Notes: This is used, for example to prevent sequential objects that are created from a parallel object; such as the KSP created by
1015     PCBJACOBI from all printing the same help messages to the screen
1016 
1017     This will not handle more complicated situations like with GASM where children may live on any subset of the parent's processes and overlap
1018 
1019 .seealso: PetscOptionsInsert(), PetscObjectSetPrintedOptions()
1020 @*/
1021 PetscErrorCode  PetscObjectInheritPrintedOptions(PetscObject pobj,PetscObject obj)
1022 {
1023   PetscErrorCode ierr;
1024   PetscMPIInt    prank,size;
1025 
1026   PetscFunctionBegin;
1027   ierr = MPI_Comm_rank(pobj->comm,&prank);CHKERRQ(ierr);
1028   ierr = MPI_Comm_size(obj->comm,&size);CHKERRQ(ierr);
1029   if (size == 1 && prank > 0) obj->optionsprinted = PETSC_TRUE;
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 /*@
1034     PetscOptionsDestroy - Destroys an option database.
1035 
1036   Input Parameter:
1037 .  options - the PetscOptions object
1038 
1039    Level: developer
1040 
1041 .seealso: PetscOptionsInsert()
1042 @*/
1043 PetscErrorCode  PetscOptionsDestroy(PetscOptions *options)
1044 {
1045   PetscErrorCode ierr;
1046 
1047   PetscFunctionBegin;
1048   ierr = PetscOptionsClear(*options);CHKERRQ(ierr);
1049   free(*options);
1050   *options = NULL;
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 PetscErrorCode  PetscOptionsDestroyDefault(void)
1055 {
1056   PetscErrorCode ierr;
1057 
1058   ierr = PetscOptionsDestroy(&defaultoptions);if (ierr) return ierr;
1059   return 0;
1060 }
1061 
1062 
1063 /*@C
1064    PetscOptionsSetValue - Sets an option name-value pair in the options
1065    database, overriding whatever is already present.
1066 
1067    Not collective, but setting values on certain processors could cause problems
1068    for parallel objects looking for options.
1069 
1070    Input Parameters:
1071 +  options - options database, use NULL for the default global database
1072 .  name - name of option, this SHOULD have the - prepended
1073 -  value - the option value (not used for all options)
1074 
1075    Level: intermediate
1076 
1077    Note:
1078    This function can be called BEFORE PetscInitialize()
1079 
1080    Only some options have values associated with them, such as
1081    -ksp_rtol tol.  Other options stand alone, such as -ksp_monitor.
1082 
1083   Developers Note: Uses malloc() directly because PETSc may not yet have been fully initialized
1084 
1085   Concepts: options database^adding option
1086 
1087 .seealso: PetscOptionsInsert()
1088 @*/
1089 PetscErrorCode  PetscOptionsSetValue(PetscOptions options,const char iname[],const char value[])
1090 {
1091   size_t         len;
1092   PetscErrorCode ierr;
1093   PetscInt       N,n,i;
1094   char           **names;
1095   char           fullname[2048];
1096   const char     *name = iname;
1097   int            match;
1098 
1099   if (!options) {
1100     if (!defaultoptions) {
1101       ierr = PetscOptionsCreateDefault();
1102       if (ierr) return ierr;
1103     }
1104     options = defaultoptions;
1105   }
1106 
1107   /* this is so that -h and -help are equivalent (p4 does not like -help)*/
1108   match = strcmp(name,"-h");
1109   if (!match) name = "-help";
1110 
1111   name++; /* skip starting hyphen */
1112   if (options->prefixind > 0) {
1113     strncpy(fullname,options->prefix,sizeof(fullname));
1114     strncat(fullname,name,sizeof(fullname)-strlen(fullname)-1);
1115     name = fullname;
1116   }
1117 
1118   /* check against aliases */
1119   N = options->Naliases;
1120   for (i=0; i<N; i++) {
1121 #if defined(PETSC_HAVE_STRCASECMP)
1122     match = strcasecmp(options->aliases1[i],name);
1123 #elif defined(PETSC_HAVE_STRICMP)
1124     match = stricmp(options->aliases1[i],name);
1125 #else
1126     Error
1127 #endif
1128     if (!match) {
1129       name = options->aliases2[i];
1130       break;
1131     }
1132   }
1133 
1134   N     = options->N;
1135   n     = N;
1136   names = options->names;
1137 
1138   for (i=0; i<N; i++) {
1139 #if defined(PETSC_HAVE_STRCASECMP)
1140     match = strcasecmp(names[i],name);
1141 #elif defined(PETSC_HAVE_STRICMP)
1142     match = stricmp(names[i],name);
1143 #else
1144     Error
1145 #endif
1146     if (!match) {
1147       if (options->values[i]) free(options->values[i]);
1148       len = value ? strlen(value) : 0;
1149       if (len) {
1150         options->values[i] = (char*)malloc((len+1)*sizeof(char));
1151         if (!options->values[i]) return PETSC_ERR_MEM;
1152         strcpy(options->values[i],value);
1153       } else options->values[i] = 0;
1154       return 0;
1155     } else if (strcmp(names[i],name) > 0) {
1156       n = i;
1157       break;
1158     }
1159   }
1160   if (N >= MAXOPTIONS) abort();
1161 
1162   /* shift remaining values down 1 */
1163   for (i=N; i>n; i--) {
1164     options->names[i]  = options->names[i-1];
1165     options->values[i] = options->values[i-1];
1166     options->used[i]   = options->used[i-1];
1167   }
1168   /* insert new name and value */
1169   len = strlen(name);
1170   options->names[n] = (char*)malloc((len+1)*sizeof(char));
1171   if (!options->names[n]) return PETSC_ERR_MEM;
1172   strcpy(options->names[n],name);
1173   len = value ? strlen(value) : 0;
1174   if (len) {
1175     options->values[n] = (char*)malloc((len+1)*sizeof(char));
1176     if (!options->values[n]) return PETSC_ERR_MEM;
1177     strcpy(options->values[n],value);
1178   } else options->values[n] = NULL;
1179   options->used[n] = PETSC_FALSE;
1180   options->N++;
1181   return 0;
1182 }
1183 
1184 /*@C
1185    PetscOptionsClearValue - Clears an option name-value pair in the options
1186    database, overriding whatever is already present.
1187 
1188    Not Collective, but setting values on certain processors could cause problems
1189    for parallel objects looking for options.
1190 
1191    Input Parameter:
1192 +  options - options database, use NULL for the default global database
1193 .  name - name of option, this SHOULD have the - prepended
1194 
1195    Level: intermediate
1196 
1197    Concepts: options database^removing option
1198 .seealso: PetscOptionsInsert()
1199 @*/
1200 PetscErrorCode  PetscOptionsClearValue(PetscOptions options,const char iname[])
1201 {
1202   PetscErrorCode ierr;
1203   PetscInt       N,n,i;
1204   char           **names,*name=(char*)iname;
1205   PetscBool      gt,match;
1206 
1207   PetscFunctionBegin;
1208   options = options ? options : defaultoptions;
1209   if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with -: Instead %s",name);
1210   name++;
1211 
1212   N     = options->N; n = 0;
1213   names = options->names;
1214 
1215   for (i=0; i<N; i++) {
1216     ierr  = PetscStrcasecmp(names[i],name,&match);CHKERRQ(ierr);
1217     ierr  = PetscStrgrt(names[i],name,&gt);CHKERRQ(ierr);
1218     if (match) {
1219       if (options->names[i])  free(options->names[i]);
1220       if (options->values[i]) free(options->values[i]);
1221       PetscOptionsMonitor(name,"");
1222       break;
1223     } else if (gt) PetscFunctionReturn(0); /* it was not listed */
1224 
1225     n++;
1226   }
1227   if (n == N) PetscFunctionReturn(0); /* it was not listed */
1228 
1229   /* shift remaining values down 1 */
1230   for (i=n; i<N-1; i++) {
1231     options->names[i]  = options->names[i+1];
1232     options->values[i] = options->values[i+1];
1233     options->used[i]   = options->used[i+1];
1234   }
1235   options->N--;
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 /*@C
1240    PetscOptionsSetAlias - Makes a key and alias for another key
1241 
1242    Not Collective, but setting values on certain processors could cause problems
1243    for parallel objects looking for options.
1244 
1245    Input Parameters:
1246 +  options - options database or NULL for default global database
1247 .  inewname - the alias
1248 -  ioldname - the name that alias will refer to
1249 
1250    Level: advanced
1251 
1252 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1253            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
1254           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1255           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1256           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1257           PetscOptionsFList(), PetscOptionsEList()
1258 @*/
1259 PetscErrorCode  PetscOptionsSetAlias(PetscOptions options,const char inewname[],const char ioldname[])
1260 {
1261   PetscErrorCode ierr;
1262   PetscInt       n = options->Naliases;
1263   size_t         len;
1264   char           *newname = (char*)inewname,*oldname = (char*)ioldname;
1265 
1266   PetscFunctionBegin;
1267   options = options ? options : defaultoptions;
1268   if (newname[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aliased must have -: Instead %s",newname);
1269   if (oldname[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aliasee must have -: Instead %s",oldname);
1270   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);
1271 
1272   newname++; oldname++;
1273   ierr = PetscStrlen(newname,&len);CHKERRQ(ierr);
1274   options->aliases1[n] = (char*)malloc((len+1)*sizeof(char));
1275   ierr = PetscStrcpy(options->aliases1[n],newname);CHKERRQ(ierr);
1276   ierr = PetscStrlen(oldname,&len);CHKERRQ(ierr);
1277   options->aliases2[n] = (char*)malloc((len+1)*sizeof(char));
1278   ierr = PetscStrcpy(options->aliases2[n],oldname);CHKERRQ(ierr);
1279   options->Naliases++;
1280   PetscFunctionReturn(0);
1281 }
1282 
1283 PetscErrorCode PetscOptionsFindPair_Private(PetscOptions options,const char pre[],const char name[],char *value[],PetscBool  *flg)
1284 {
1285   PetscErrorCode ierr;
1286   PetscInt       i,N;
1287   size_t         len;
1288   char           **names,tmp[256];
1289   PetscBool      match;
1290 
1291   PetscFunctionBegin;
1292   options = options ? options : defaultoptions;
1293   N     = options->N;
1294   names = options->names;
1295 
1296   if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with -: Instead %s",name);
1297 
1298   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1299   if (pre) {
1300     char       *ptr   = tmp;
1301     const char *namep = name;
1302     if (pre[0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix should not begin with a -");
1303     if (name[1] == '-') {
1304       *ptr++ = '-';
1305       namep++;
1306     }
1307     ierr = PetscStrncpy(ptr,pre,tmp+sizeof(tmp)-ptr);CHKERRQ(ierr);
1308     tmp[sizeof(tmp)-1] = 0;
1309     ierr = PetscStrlen(tmp,&len);CHKERRQ(ierr);
1310     ierr = PetscStrncat(tmp,namep+1,sizeof(tmp)-len-1);CHKERRQ(ierr);
1311   } else {
1312     ierr = PetscStrncpy(tmp,name+1,sizeof(tmp));CHKERRQ(ierr);
1313     tmp[sizeof(tmp)-1] = 0;
1314   }
1315 #if defined(PETSC_USE_DEBUG)
1316   {
1317     PetscBool valid;
1318     char      key[sizeof(tmp)+1] = "-";
1319 
1320     ierr = PetscMemcpy(key+1,tmp,sizeof(tmp));CHKERRQ(ierr);
1321     ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr);
1322     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1323   }
1324 #endif
1325 
1326   /* slow search */
1327   *flg = PETSC_FALSE;
1328   for (i=0; i<N; i++) {
1329     ierr = PetscStrcasecmp(names[i],tmp,&match);CHKERRQ(ierr);
1330     if (match) {
1331       *value           = options->values[i];
1332       options->used[i] = PETSC_TRUE;
1333       *flg             = PETSC_TRUE;
1334       break;
1335     }
1336   }
1337   if (!*flg) {
1338     PetscInt j,cnt = 0,locs[16],loce[16];
1339     size_t   n;
1340     ierr = PetscStrlen(tmp,&n);CHKERRQ(ierr);
1341     /* determine the location and number of all _%d_ in the key */
1342     for (i=0; i< (PetscInt)n; i++) {
1343       if (tmp[i] == '_') {
1344         for (j=i+1; j< (PetscInt)n; j++) {
1345           if (tmp[j] >= '0' && tmp[j] <= '9') continue;
1346           if (tmp[j] == '_' && j > i+1) { /* found a number */
1347             locs[cnt]   = i+1;
1348             loce[cnt++] = j+1;
1349           }
1350           break;
1351         }
1352       }
1353     }
1354     if (cnt) {
1355       char tmp2[256];
1356       for (i=0; i<cnt; i++) {
1357         ierr = PetscStrcpy(tmp2,"-");CHKERRQ(ierr);
1358         ierr = PetscStrncat(tmp2,tmp,locs[i]);CHKERRQ(ierr);
1359         ierr = PetscStrcat(tmp2,tmp+loce[i]);CHKERRQ(ierr);
1360         ierr = PetscOptionsFindPair_Private(options,NULL,tmp2,value,flg);CHKERRQ(ierr);
1361         if (*flg) break;
1362       }
1363     }
1364   }
1365   PetscFunctionReturn(0);
1366 }
1367 
1368 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[], char *value[], PetscBool *flg)
1369 {
1370   PetscErrorCode ierr;
1371   PetscInt       i,N;
1372   size_t         len;
1373   char           **names,tmp[256];
1374   PetscBool      match;
1375 
1376   PetscFunctionBegin;
1377   options = options ? options : defaultoptions;
1378   N     = options->N;
1379   names = options->names;
1380 
1381   if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with -: Instead %s",name);
1382 
1383   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1384   if (pre) {
1385     char       *ptr   = tmp;
1386     const char *namep = name;
1387     if (pre[0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix should not begin with a -");
1388     if (name[1] == '-') {
1389       *ptr++ = '-';
1390       namep++;
1391     }
1392     ierr = PetscStrncpy(ptr,pre,tmp+sizeof(tmp)-ptr);CHKERRQ(ierr);
1393     tmp[sizeof(tmp)-1] = 0;
1394     ierr = PetscStrlen(tmp,&len);CHKERRQ(ierr);
1395     ierr = PetscStrncat(tmp,namep+1,sizeof(tmp)-len-1);CHKERRQ(ierr);
1396   } else {
1397     ierr = PetscStrncpy(tmp,name+1,sizeof(tmp));CHKERRQ(ierr);
1398     tmp[sizeof(tmp)-1] = 0;
1399   }
1400 #if defined(PETSC_USE_DEBUG)
1401   {
1402     PetscBool valid;
1403     char      key[sizeof(tmp)+1] = "-";
1404 
1405     ierr = PetscMemcpy(key+1,tmp,sizeof(tmp));CHKERRQ(ierr);
1406     ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr);
1407     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1408   }
1409 #endif
1410 
1411   /* slow search */
1412   *flg = PETSC_FALSE;
1413   ierr = PetscStrlen(tmp,&len);CHKERRQ(ierr);
1414   for (i = 0; i < N; ++i) {
1415     ierr = PetscStrncmp(names[i], tmp, len, &match);CHKERRQ(ierr);
1416     if (match) {
1417       if (value) *value = options->values[i];
1418       options->used[i]  = PETSC_TRUE;
1419       if (flg)   *flg   = PETSC_TRUE;
1420       break;
1421     }
1422   }
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 /*@C
1427    PetscOptionsReject - Generates an error if a certain option is given.
1428 
1429    Not Collective, but setting values on certain processors could cause problems
1430    for parallel objects looking for options.
1431 
1432    Input Parameters:
1433 +  options - options database use NULL for default global database
1434 .  name - the option one is seeking
1435 -  mess - error message (may be NULL)
1436 
1437    Level: advanced
1438 
1439    Concepts: options database^rejecting option
1440 
1441 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1442            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1443           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1444           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1445           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1446           PetscOptionsFList(), PetscOptionsEList()
1447 @*/
1448 PetscErrorCode  PetscOptionsReject(PetscOptions options,const char name[],const char mess[])
1449 {
1450   PetscErrorCode ierr;
1451   PetscBool      flag = PETSC_FALSE;
1452 
1453   PetscFunctionBegin;
1454   options = options ? options : defaultoptions;
1455   ierr = PetscOptionsHasName(options,NULL,name,&flag);CHKERRQ(ierr);
1456   if (flag) {
1457     if (mess) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: %s with %s",name,mess);
1458     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: %s",name);
1459   }
1460   PetscFunctionReturn(0);
1461 }
1462 
1463 /*@C
1464    PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or boolean, even
1465                       its value is set to false.
1466 
1467    Not Collective
1468 
1469    Input Parameters:
1470 +  options - options database use NULL for default global database
1471 .  name - the option one is seeking
1472 -  pre - string to prepend to the name or NULL
1473 
1474    Output Parameters:
1475 .  set - PETSC_TRUE if found else PETSC_FALSE.
1476 
1477    Level: beginner
1478 
1479    Concepts: options database^has option name
1480 
1481    Notes: Name cannot be simply -h
1482 
1483           In many cases you probably want to use PetscOptionsGetBool() instead of calling this, to allowing toggling values.
1484 
1485 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1486            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1487           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1488           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1489           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1490           PetscOptionsFList(), PetscOptionsEList()
1491 @*/
1492 PetscErrorCode  PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool  *set)
1493 {
1494   char           *value;
1495   PetscErrorCode ierr;
1496   PetscBool      flag;
1497 
1498   PetscFunctionBegin;
1499   options = options ? options : defaultoptions;
1500   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1501   if (set) *set = flag;
1502   PetscFunctionReturn(0);
1503 }
1504 
1505 /*@C
1506    PetscOptionsGetInt - Gets the integer value for a particular option in the database.
1507 
1508    Not Collective
1509 
1510    Input Parameters:
1511 +  options - options database use NULL for default global database
1512 .  pre - the string to prepend to the name or NULL
1513 -  name - the option one is seeking
1514 
1515    Output Parameter:
1516 +  ivalue - the integer value to return
1517 -  set - PETSC_TRUE if found, else PETSC_FALSE
1518 
1519    Level: beginner
1520 
1521    Concepts: options database^has int
1522 
1523 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1524           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1525           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1526           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1527           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1528           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1529           PetscOptionsFList(), PetscOptionsEList()
1530 @*/
1531 PetscErrorCode  PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool  *set)
1532 {
1533   char           *value;
1534   PetscErrorCode ierr;
1535   PetscBool      flag;
1536 
1537   PetscFunctionBegin;
1538   PetscValidCharPointer(name,2);
1539   PetscValidIntPointer(ivalue,3);
1540   options = options ? options : defaultoptions;
1541   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1542   if (flag) {
1543     if (!value) {
1544       if (set) *set = PETSC_FALSE;
1545     } else {
1546       if (set) *set = PETSC_TRUE;
1547       ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr);
1548     }
1549   } else {
1550     if (set) *set = PETSC_FALSE;
1551   }
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 /*@C
1556      PetscOptionsGetEList - Puts a list of option values that a single one may be selected from
1557 
1558    Not Collective
1559 
1560    Input Parameters:
1561 +  options - options database use NULL for default global database
1562 .  pre - the string to prepend to the name or NULL
1563 .  opt - option name
1564 .  list - the possible choices (one of these must be selected, anything else is invalid)
1565 .  ntext - number of choices
1566 
1567    Output Parameter:
1568 +  value - the index of the value to return (defaults to zero if the option name is given but choice is listed)
1569 -  set - PETSC_TRUE if found, else PETSC_FALSE
1570 
1571    Level: intermediate
1572 
1573    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
1574 
1575    Concepts: options database^list
1576 
1577 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1578            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1579           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1580           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1581           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1582           PetscOptionsFList(), PetscOptionsEList()
1583 @*/
1584 PetscErrorCode  PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool  *set)
1585 {
1586   PetscErrorCode ierr;
1587   size_t         alen,len = 0;
1588   char           *svalue;
1589   PetscBool      aset,flg = PETSC_FALSE;
1590   PetscInt       i;
1591 
1592   PetscFunctionBegin;
1593   options = options ? options : defaultoptions;
1594   for (i=0; i<ntext; i++) {
1595     ierr = PetscStrlen(list[i],&alen);CHKERRQ(ierr);
1596     if (alen > len) len = alen;
1597   }
1598   len += 5; /* a little extra space for user mistypes */
1599   ierr = PetscMalloc1(len,&svalue);CHKERRQ(ierr);
1600   ierr = PetscOptionsGetString(options,pre,opt,svalue,len,&aset);CHKERRQ(ierr);
1601   if (aset) {
1602     ierr = PetscEListFind(ntext,list,svalue,value,&flg);CHKERRQ(ierr);
1603     if (!flg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s",svalue,pre ? pre : "",opt+1);
1604     if (set) *set = PETSC_TRUE;
1605   } else if (set) *set = PETSC_FALSE;
1606   ierr = PetscFree(svalue);CHKERRQ(ierr);
1607   PetscFunctionReturn(0);
1608 }
1609 
1610 /*@C
1611    PetscOptionsGetEnum - Gets the enum value for a particular option in the database.
1612 
1613    Not Collective
1614 
1615    Input Parameters:
1616 +  options - options database use NULL for default global database
1617 .  pre - option prefix or NULL
1618 .  opt - option name
1619 .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
1620 -  defaultv - the default (current) value
1621 
1622    Output Parameter:
1623 +  value - the  value to return
1624 -  set - PETSC_TRUE if found, else PETSC_FALSE
1625 
1626    Level: beginner
1627 
1628    Concepts: options database
1629 
1630    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1631 
1632           list is usually something like PCASMTypes or some other predefined list of enum names
1633 
1634 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1635           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1636           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1637           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1638           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1639           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1640           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
1641 @*/
1642 PetscErrorCode  PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool  *set)
1643 {
1644   PetscErrorCode ierr;
1645   PetscInt       ntext = 0,tval;
1646   PetscBool      fset;
1647 
1648   PetscFunctionBegin;
1649   options = options ? options : defaultoptions;
1650   while (list[ntext++]) {
1651     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
1652   }
1653   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
1654   ntext -= 3;
1655   ierr = PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);CHKERRQ(ierr);
1656   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
1657   if (fset) *value = (PetscEnum)tval;
1658   if (set) *set = fset;
1659   PetscFunctionReturn(0);
1660 }
1661 
1662 /*@C
1663    PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
1664             option in the database.
1665 
1666    Not Collective
1667 
1668    Input Parameters:
1669 +  options - options database use NULL for default global database
1670 .  pre - the string to prepend to the name or NULL
1671 -  name - the option one is seeking
1672 
1673    Output Parameter:
1674 +  ivalue - the logical value to return
1675 -  set - PETSC_TRUE  if found, else PETSC_FALSE
1676 
1677    Level: beginner
1678 
1679    Notes:
1680        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
1681        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
1682 
1683        If the user does not supply the option (as either true or false) ivalue is NOT changed. Thus
1684      you NEED TO ALWAYS initialize the ivalue.
1685 
1686    Concepts: options database^has logical
1687 
1688 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1689           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(),
1690           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1691           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1692           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1693           PetscOptionsFList(), PetscOptionsEList()
1694 @*/
1695 PetscErrorCode  PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool  *ivalue,PetscBool  *set)
1696 {
1697   char           *value;
1698   PetscBool      flag;
1699   PetscErrorCode ierr;
1700 
1701   PetscFunctionBegin;
1702   PetscValidCharPointer(name,2);
1703   if (ivalue) PetscValidIntPointer(ivalue,3);
1704   options = options ? options : defaultoptions;
1705   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1706   if (flag) {
1707     if (set) *set = PETSC_TRUE;
1708     if (!value) {
1709       if (ivalue) *ivalue = PETSC_TRUE;
1710     } else {
1711       ierr = PetscOptionsStringToBool(value, ivalue);CHKERRQ(ierr);
1712     }
1713   } else {
1714     if (set) *set = PETSC_FALSE;
1715   }
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 /*@C
1720    PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
1721    option in the database.  The values must be separated with commas with
1722    no intervening spaces.
1723 
1724    Not Collective
1725 
1726    Input Parameters:
1727 +  options - options database use NULL for default global database
1728 .  pre - string to prepend to each name or NULL
1729 .  name - the option one is seeking
1730 -  nmax - maximum number of values to retrieve
1731 
1732    Output Parameter:
1733 +  dvalue - the integer values to return
1734 .  nmax - actual number of values retreived
1735 -  set - PETSC_TRUE if found, else PETSC_FALSE
1736 
1737    Level: beginner
1738 
1739    Concepts: options database^array of ints
1740 
1741    Notes:
1742        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
1743        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
1744 
1745 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
1746            PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1747           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1748           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1749           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1750           PetscOptionsFList(), PetscOptionsEList()
1751 @*/
1752 PetscErrorCode  PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool  *set)
1753 {
1754   char           *value;
1755   PetscErrorCode ierr;
1756   PetscInt       n = 0;
1757   PetscBool      flag;
1758   PetscToken     token;
1759 
1760   PetscFunctionBegin;
1761   PetscValidCharPointer(name,2);
1762   PetscValidIntPointer(dvalue,3);
1763   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1764   if (!flag)  {if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);}
1765   if (!value) {if (set) *set = PETSC_TRUE;  *nmax = 0; PetscFunctionReturn(0);}
1766 
1767   if (set) *set = PETSC_TRUE;
1768 
1769   ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
1770   ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
1771   while (n < *nmax) {
1772     if (!value) break;
1773     ierr = PetscOptionsStringToBool(value,dvalue);CHKERRQ(ierr);
1774     ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
1775     dvalue++;
1776     n++;
1777   }
1778   ierr  = PetscTokenDestroy(&token);CHKERRQ(ierr);
1779   *nmax = n;
1780   PetscFunctionReturn(0);
1781 }
1782 
1783 /*@C
1784    PetscOptionsGetReal - Gets the double precision value for a particular
1785    option in the database.
1786 
1787    Not Collective
1788 
1789    Input Parameters:
1790 +  options - options database use NULL for default global database
1791 .  pre - string to prepend to each name or NULL
1792 -  name - the option one is seeking
1793 
1794    Output Parameter:
1795 +  dvalue - the double value to return
1796 -  set - PETSC_TRUE if found, PETSC_FALSE if not found
1797 
1798    Note: if the option is given but no value is provided then set is given the value PETSC_FALSE
1799 
1800    Level: beginner
1801 
1802    Concepts: options database^has double
1803 
1804 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
1805            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
1806           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1807           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1808           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1809           PetscOptionsFList(), PetscOptionsEList()
1810 @*/
1811 PetscErrorCode  PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool  *set)
1812 {
1813   char           *value;
1814   PetscErrorCode ierr;
1815   PetscBool      flag;
1816 
1817   PetscFunctionBegin;
1818   PetscValidCharPointer(name,2);
1819   PetscValidRealPointer(dvalue,3);
1820   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1821   if (flag) {
1822     if (!value) {
1823       if (set) *set = PETSC_FALSE;
1824     } else {
1825       if (set) *set = PETSC_TRUE;
1826       ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
1827     }
1828   } else {
1829     if (set) *set = PETSC_FALSE;
1830   }
1831   PetscFunctionReturn(0);
1832 }
1833 
1834 /*@C
1835    PetscOptionsGetScalar - Gets the scalar value for a particular
1836    option in the database.
1837 
1838    Not Collective
1839 
1840    Input Parameters:
1841 +  options - options database use NULL for default global database
1842 .  pre - string to prepend to each name or NULL
1843 -  name - the option one is seeking
1844 
1845    Output Parameter:
1846 +  dvalue - the double value to return
1847 -  set - PETSC_TRUE if found, else PETSC_FALSE
1848 
1849    Level: beginner
1850 
1851    Usage:
1852    A complex number 2+3i must be specified with NO spaces
1853 
1854    Note: if the option is given but no value is provided then set is given the value PETSC_FALSE
1855 
1856    Concepts: options database^has scalar
1857 
1858 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
1859            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1860           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1861           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1862           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1863           PetscOptionsFList(), PetscOptionsEList()
1864 @*/
1865 PetscErrorCode  PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool  *set)
1866 {
1867   char           *value;
1868   PetscBool      flag;
1869   PetscErrorCode ierr;
1870 
1871   PetscFunctionBegin;
1872   PetscValidCharPointer(name,2);
1873   PetscValidScalarPointer(dvalue,3);
1874   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1875   if (flag) {
1876     if (!value) {
1877       if (set) *set = PETSC_FALSE;
1878     } else {
1879 #if !defined(PETSC_USE_COMPLEX)
1880       ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
1881 #else
1882       ierr = PetscOptionsStringToScalar(value,dvalue);CHKERRQ(ierr);
1883 #endif
1884       if (set) *set = PETSC_TRUE;
1885     }
1886   } else { /* flag */
1887     if (set) *set = PETSC_FALSE;
1888   }
1889   PetscFunctionReturn(0);
1890 }
1891 
1892 /*@C
1893    PetscOptionsGetRealArray - Gets an array of double precision values for a
1894    particular option in the database.  The values must be separated with
1895    commas with no intervening spaces.
1896 
1897    Not Collective
1898 
1899    Input Parameters:
1900 +  options - options database use NULL for default global database
1901 .  pre - string to prepend to each name or NULL
1902 .  name - the option one is seeking
1903 -  nmax - maximum number of values to retrieve
1904 
1905    Output Parameters:
1906 +  dvalue - the double values to return
1907 .  nmax - actual number of values retreived
1908 -  set - PETSC_TRUE if found, else PETSC_FALSE
1909 
1910    Level: beginner
1911 
1912    Concepts: options database^array of doubles
1913 
1914 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
1915            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
1916           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1917           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1918           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1919           PetscOptionsFList(), PetscOptionsEList()
1920 @*/
1921 PetscErrorCode  PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool  *set)
1922 {
1923   char           *value;
1924   PetscErrorCode ierr;
1925   PetscInt       n = 0;
1926   PetscBool      flag;
1927   PetscToken     token;
1928 
1929   PetscFunctionBegin;
1930   PetscValidCharPointer(name,2);
1931   PetscValidRealPointer(dvalue,3);
1932   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
1933   if (!flag) {
1934     if (set) *set = PETSC_FALSE;
1935     *nmax = 0;
1936     PetscFunctionReturn(0);
1937   }
1938   if (!value) {
1939     if (set) *set = PETSC_TRUE;
1940     *nmax = 0;
1941     PetscFunctionReturn(0);
1942   }
1943 
1944   if (set) *set = PETSC_TRUE;
1945 
1946   ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
1947   ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
1948   while (n < *nmax) {
1949     if (!value) break;
1950     ierr = PetscOptionsStringToReal(value,dvalue++);CHKERRQ(ierr);
1951     ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
1952     n++;
1953   }
1954   ierr  = PetscTokenDestroy(&token);CHKERRQ(ierr);
1955   *nmax = n;
1956   PetscFunctionReturn(0);
1957 }
1958 
1959 /*@C
1960    PetscOptionsGetScalarArray - Gets an array of scalars for a
1961    particular option in the database.  The values must be separated with
1962    commas with no intervening spaces.
1963 
1964    Not Collective
1965 
1966    Input Parameters:
1967 +  options - options database use NULL for default global database
1968 .  pre - string to prepend to each name or NULL
1969 .  name - the option one is seeking
1970 -  nmax - maximum number of values to retrieve
1971 
1972    Output Parameters:
1973 +  dvalue - the scalar values to return
1974 .  nmax - actual number of values retreived
1975 -  set - PETSC_TRUE if found, else PETSC_FALSE
1976 
1977    Level: beginner
1978 
1979    Concepts: options database^array of doubles
1980 
1981 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
1982            PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
1983           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1984           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1985           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1986           PetscOptionsFList(), PetscOptionsEList()
1987 @*/
1988 PetscErrorCode  PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],PetscInt *nmax,PetscBool  *set)
1989 {
1990   char           *value;
1991   PetscErrorCode ierr;
1992   PetscInt       n = 0;
1993   PetscBool      flag;
1994   PetscToken     token;
1995 
1996   PetscFunctionBegin;
1997   PetscValidCharPointer(name,2);
1998   PetscValidRealPointer(dvalue,3);
1999   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
2000   if (!flag) {
2001     if (set) *set = PETSC_FALSE;
2002     *nmax = 0;
2003     PetscFunctionReturn(0);
2004   }
2005   if (!value) {
2006     if (set) *set = PETSC_TRUE;
2007     *nmax = 0;
2008     PetscFunctionReturn(0);
2009   }
2010 
2011   if (set) *set = PETSC_TRUE;
2012 
2013   ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
2014   ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2015   while (n < *nmax) {
2016     if (!value) break;
2017     ierr = PetscOptionsStringToScalar(value,dvalue++);CHKERRQ(ierr);
2018     ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2019     n++;
2020   }
2021   ierr  = PetscTokenDestroy(&token);CHKERRQ(ierr);
2022   *nmax = n;
2023   PetscFunctionReturn(0);
2024 }
2025 
2026 /*@C
2027    PetscOptionsGetIntArray - Gets an array of integer values for a particular
2028    option in the database.
2029 
2030    Not Collective
2031 
2032    Input Parameters:
2033 +  options - options database use NULL for default global database
2034 .  pre - string to prepend to each name or NULL
2035 .  name - the option one is seeking
2036 -  nmax - maximum number of values to retrieve
2037 
2038    Output Parameter:
2039 +  dvalue - the integer values to return
2040 .  nmax - actual number of values retreived
2041 -  set - PETSC_TRUE if found, else PETSC_FALSE
2042 
2043    Level: beginner
2044 
2045    Notes:
2046    The array can be passed as
2047    a comma separated list:                                 0,1,2,3,4,5,6,7
2048    a range (start-end+1):                                  0-8
2049    a range with given increment (start-end+1:inc):         0-7:2
2050    a combination of values and ranges separated by commas: 0,1-8,8-15:2
2051 
2052    There must be no intervening spaces between the values.
2053 
2054    Concepts: options database^array of ints
2055 
2056 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2057            PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2058           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2059           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2060           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2061           PetscOptionsFList(), PetscOptionsEList()
2062 @*/
2063 PetscErrorCode  PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt dvalue[],PetscInt *nmax,PetscBool  *set)
2064 {
2065   char           *value;
2066   PetscErrorCode ierr;
2067   PetscInt       n = 0,i,j,start,end,inc,nvalues;
2068   size_t         len;
2069   PetscBool      flag,foundrange;
2070   PetscToken     token;
2071 
2072   PetscFunctionBegin;
2073   PetscValidCharPointer(name,2);
2074   PetscValidIntPointer(dvalue,3);
2075   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
2076   if (!flag) {
2077     if (set) *set = PETSC_FALSE;
2078     *nmax = 0;
2079     PetscFunctionReturn(0);
2080   }
2081   if (!value) {
2082     if (set) *set = PETSC_TRUE;
2083     *nmax = 0;
2084     PetscFunctionReturn(0);
2085   }
2086 
2087   if (set) *set = PETSC_TRUE;
2088 
2089   ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
2090   ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2091   while (n < *nmax) {
2092     if (!value) break;
2093 
2094     /* look for form  d-D where d and D are integers */
2095     foundrange = PETSC_FALSE;
2096     ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
2097     if (value[0] == '-') i=2;
2098     else i=1;
2099     for (;i<(int)len; i++) {
2100       if (value[i] == '-') {
2101         if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
2102         value[i] = 0;
2103 
2104         ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
2105         inc  = 1;
2106         j    = i+1;
2107         for (;j<(int)len; j++) {
2108           if (value[j] == ':') {
2109             value[j] = 0;
2110 
2111             ierr = PetscOptionsStringToInt(value+j+1,&inc);CHKERRQ(ierr);
2112             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);CHKERRQ(ierr);
2113             break;
2114           }
2115         }
2116         ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
2117         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);
2118         nvalues = (end-start)/inc + (end-start)%inc;
2119         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);
2120         for (;start<end; start+=inc) {
2121           *dvalue = start; dvalue++;n++;
2122         }
2123         foundrange = PETSC_TRUE;
2124         break;
2125       }
2126     }
2127     if (!foundrange) {
2128       ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
2129       dvalue++;
2130       n++;
2131     }
2132     ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2133   }
2134   ierr  = PetscTokenDestroy(&token);CHKERRQ(ierr);
2135   *nmax = n;
2136   PetscFunctionReturn(0);
2137 }
2138 
2139 /*@C
2140    PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.
2141 
2142    Not Collective
2143 
2144    Input Parameters:
2145 +  options - options database use NULL for default global database
2146 .  pre - option prefix or NULL
2147 .  name - option name
2148 .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2149 -  nmax - maximum number of values to retrieve
2150 
2151    Output Parameters:
2152 +  dvalue - the  enum values to return
2153 .  nmax - actual number of values retreived
2154 -  set - PETSC_TRUE if found, else PETSC_FALSE
2155 
2156    Level: beginner
2157 
2158    Concepts: options database
2159 
2160    Notes:
2161    The array must be passed as a comma separated list.
2162 
2163    There must be no intervening spaces between the values.
2164 
2165    list is usually something like PCASMTypes or some other predefined list of enum names.
2166 
2167 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2168           PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2169           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(),
2170           PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(),
2171           PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2172           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2173 @*/
2174 PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum dvalue[],PetscInt *nmax,PetscBool *set)
2175 {
2176   char           *svalue;
2177   PetscInt       n = 0;
2178   PetscEnum      evalue;
2179   PetscBool      flag;
2180   PetscToken     token;
2181   PetscErrorCode ierr;
2182 
2183   PetscFunctionBegin;
2184   PetscValidCharPointer(name,2);
2185   PetscValidPointer(list,3);
2186   PetscValidPointer(dvalue,4);
2187   PetscValidIntPointer(nmax,5);
2188 
2189   ierr = PetscOptionsFindPair_Private(options,pre,name,&svalue,&flag);CHKERRQ(ierr);
2190   if (!flag) {
2191     if (set) *set = PETSC_FALSE;
2192     *nmax = 0;
2193     PetscFunctionReturn(0);
2194   }
2195   if (!svalue) {
2196     if (set) *set = PETSC_TRUE;
2197     *nmax = 0;
2198     PetscFunctionReturn(0);
2199   }
2200   if (set) *set = PETSC_TRUE;
2201 
2202   ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr);
2203   ierr = PetscTokenFind(token,&svalue);CHKERRQ(ierr);
2204   while (svalue && n < *nmax) {
2205     ierr = PetscEnumFind(list,svalue,&evalue,&flag);CHKERRQ(ierr);
2206     if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1);
2207     dvalue[n++] = evalue;
2208     ierr = PetscTokenFind(token,&svalue);CHKERRQ(ierr);
2209   }
2210   *nmax = n;
2211   ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
2212   PetscFunctionReturn(0);
2213 }
2214 
2215 /*@C
2216    PetscOptionsGetString - Gets the string value for a particular option in
2217    the database.
2218 
2219    Not Collective
2220 
2221    Input Parameters:
2222 +  options - options database use NULL for default global database
2223 .  pre - string to prepend to name or NULL
2224 .  name - the option one is seeking
2225 -  len - maximum length of the string including null termination
2226 
2227    Output Parameters:
2228 +  string - location to copy string
2229 -  set - PETSC_TRUE if found, else PETSC_FALSE
2230 
2231    Level: beginner
2232 
2233    Fortran Note:
2234    The Fortran interface is slightly different from the C/C++
2235    interface (len is not used).  Sample usage in Fortran follows
2236 .vb
2237       character *20    string
2238       PetscErrorCode   ierr
2239       PetscBool        set
2240       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2241 .ve
2242 
2243    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
2244 
2245    Concepts: options database^string
2246 
2247     Note:
2248       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).
2249 
2250 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2251            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2252           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2253           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2254           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2255           PetscOptionsFList(), PetscOptionsEList()
2256 @*/
2257 PetscErrorCode  PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool  *set)
2258 {
2259   char           *value;
2260   PetscErrorCode ierr;
2261   PetscBool      flag;
2262 
2263   PetscFunctionBegin;
2264   PetscValidCharPointer(name,2);
2265   PetscValidCharPointer(string,3);
2266   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
2267   if (!flag) {
2268     if (set) *set = PETSC_FALSE;
2269   } else {
2270     if (set) *set = PETSC_TRUE;
2271     if (value) {
2272       ierr = PetscStrncpy(string,value,len);CHKERRQ(ierr);
2273       string[len-1] = 0;        /* Ensure that the string is NULL terminated */
2274     } else {
2275       ierr = PetscMemzero(string,len);CHKERRQ(ierr);
2276     }
2277   }
2278   PetscFunctionReturn(0);
2279 }
2280 
2281 char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[])
2282 {
2283   char           *value;
2284   PetscErrorCode ierr;
2285   PetscBool      flag;
2286 
2287   PetscFunctionBegin;
2288   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(0);
2289   if (flag) PetscFunctionReturn(value);
2290   else PetscFunctionReturn(0);
2291 }
2292 
2293 
2294 /*@C
2295    PetscOptionsGetStringArray - Gets an array of string values for a particular
2296    option in the database. The values must be separated with commas with
2297    no intervening spaces.
2298 
2299    Not Collective
2300 
2301    Input Parameters:
2302 +  options - options database use NULL for default global database
2303 .  pre - string to prepend to name or NULL
2304 .  name - the option one is seeking
2305 -  nmax - maximum number of strings
2306 
2307    Output Parameter:
2308 +  strings - location to copy strings
2309 -  set - PETSC_TRUE if found, else PETSC_FALSE
2310 
2311    Level: beginner
2312 
2313    Notes:
2314    The user should pass in an array of pointers to char, to hold all the
2315    strings returned by this function.
2316 
2317    The user is responsible for deallocating the strings that are
2318    returned. The Fortran interface for this routine is not supported.
2319 
2320    Contributed by Matthew Knepley.
2321 
2322    Concepts: options database^array of strings
2323 
2324 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2325            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2326           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2327           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2328           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2329           PetscOptionsFList(), PetscOptionsEList()
2330 @*/
2331 PetscErrorCode  PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool  *set)
2332 {
2333   char           *value;
2334   PetscErrorCode ierr;
2335   PetscInt       n;
2336   PetscBool      flag;
2337   PetscToken     token;
2338 
2339   PetscFunctionBegin;
2340   PetscValidCharPointer(name,2);
2341   PetscValidPointer(strings,3);
2342   ierr = PetscOptionsFindPair_Private(options,pre,name,&value,&flag);CHKERRQ(ierr);
2343   if (!flag) {
2344     *nmax = 0;
2345     if (set) *set = PETSC_FALSE;
2346     PetscFunctionReturn(0);
2347   }
2348   if (!value) {
2349     *nmax = 0;
2350     if (set) *set = PETSC_FALSE;
2351     PetscFunctionReturn(0);
2352   }
2353   if (!*nmax) {
2354     if (set) *set = PETSC_FALSE;
2355     PetscFunctionReturn(0);
2356   }
2357   if (set) *set = PETSC_TRUE;
2358 
2359   ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
2360   ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2361   n    = 0;
2362   while (n < *nmax) {
2363     if (!value) break;
2364     ierr = PetscStrallocpy(value,&strings[n]);CHKERRQ(ierr);
2365     ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
2366     n++;
2367   }
2368   ierr  = PetscTokenDestroy(&token);CHKERRQ(ierr);
2369   *nmax = n;
2370   PetscFunctionReturn(0);
2371 }
2372 
2373 /*@C
2374    PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database
2375 
2376    Not Collective
2377 
2378    Input Parameter:
2379 +   options - options database use NULL for default global database
2380 -   option - string name of option
2381 
2382    Output Parameter:
2383 .   used - PETSC_TRUE if the option was used, otherwise false, including if option was not found in options database
2384 
2385    Level: advanced
2386 
2387 .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed()
2388 @*/
2389 PetscErrorCode  PetscOptionsUsed(PetscOptions options,const char *option,PetscBool *used)
2390 {
2391   PetscInt       i;
2392   PetscErrorCode ierr;
2393 
2394   PetscFunctionBegin;
2395   options = options ? options : defaultoptions;
2396   *used = PETSC_FALSE;
2397   for (i=0; i<options->N; i++) {
2398     ierr = PetscStrcmp(options->names[i],option,used);CHKERRQ(ierr);
2399     if (*used) {
2400       *used = options->used[i];
2401       break;
2402     }
2403   }
2404   PetscFunctionReturn(0);
2405 }
2406 
2407 /*@C
2408    PetscOptionsAllUsed - Returns a count of the number of options in the
2409    database that have never been selected.
2410 
2411    Not Collective
2412 
2413    Input Parameter:
2414 .  options - options database use NULL for default global database
2415 
2416    Output Parameter:
2417 .   N - count of options not used
2418 
2419    Level: advanced
2420 
2421 .seealso: PetscOptionsView()
2422 @*/
2423 PetscErrorCode  PetscOptionsAllUsed(PetscOptions options,PetscInt *N)
2424 {
2425   PetscInt     i,n = 0;
2426 
2427   PetscFunctionBegin;
2428   options = options ? options : defaultoptions;
2429   for (i=0; i<options->N; i++) {
2430     if (!options->used[i]) n++;
2431   }
2432   *N = n;
2433   PetscFunctionReturn(0);
2434 }
2435 
2436 /*@C
2437     PetscOptionsLeft - Prints to screen any options that were set and never used.
2438 
2439   Not collective
2440 
2441    Input Parameter:
2442 .  options - options database use NULL for default global database
2443 
2444    Options Database Key:
2445 .  -options_left - Activates OptionsAllUsed() within PetscFinalize()
2446 
2447   Level: advanced
2448 
2449 .seealso: PetscOptionsAllUsed()
2450 @*/
2451 PetscErrorCode  PetscOptionsLeft(PetscOptions options)
2452 {
2453   PetscErrorCode ierr;
2454   PetscInt       i;
2455 
2456   PetscFunctionBegin;
2457   options = options ? options : defaultoptions;
2458   for (i=0; i<options->N; i++) {
2459     if (!options->used[i]) {
2460       if (options->values[i]) {
2461         ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",options->names[i],options->values[i]);CHKERRQ(ierr);
2462       } else {
2463         ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",options->names[i]);CHKERRQ(ierr);
2464       }
2465     }
2466   }
2467   PetscFunctionReturn(0);
2468 }
2469 
2470 /*@C
2471   PetscOptionsLeftGet - Returns all options that were set and never used.
2472 
2473   Not collective
2474 
2475    Input Parameter:
2476 .  options - options database use NULL for default global database
2477 
2478    Output Parameter:
2479 .   N - count of options not used
2480 .   names - names of options not used
2481 .   values - values of options not used
2482 
2483   Level: advanced
2484 
2485   Notes:
2486   Users should call PetscOptionsLeftRestore() to free the memory allocated in this routine
2487 
2488 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft()
2489 @*/
2490 PetscErrorCode  PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[])
2491 {
2492   PetscErrorCode ierr;
2493   PetscInt       i,n;
2494 
2495   PetscFunctionBegin;
2496   options = options ? options : defaultoptions;
2497 
2498   /* The number of unused PETSc options */
2499   n = 0;
2500   for (i=0; i<options->N; i++) {
2501     if (!options->used[i]) {
2502       n++;
2503     }
2504   }
2505   if (N) {*N = n;}
2506   if (names)  { ierr = PetscMalloc1(n,names);CHKERRQ(ierr); }
2507   if (values) { ierr = PetscMalloc1(n,values);CHKERRQ(ierr); }
2508 
2509   n = 0;
2510   if (names || values) {
2511     for (i=0; i<options->N; i++) {
2512       if (!options->used[i]) {
2513         if (names)  (*names)[n]  = options->names[i];
2514         if (values) (*values)[n] = options->values[i];
2515         n++;
2516       }
2517     }
2518   }
2519   PetscFunctionReturn(0);
2520 }
2521 
2522 
2523 /*@C
2524   PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using PetscOptionsLeftGet.
2525 
2526   Not collective
2527 
2528    Input Parameter:
2529 .   options - options database use NULL for default global database
2530 .   names - names of options not used
2531 .   values - values of options not used
2532 
2533   Level: advanced
2534 
2535 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet
2536 @*/
2537 PetscErrorCode  PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[])
2538 {
2539   PetscErrorCode ierr;
2540 
2541   PetscFunctionBegin;
2542   options = options ? options : defaultoptions;
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