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