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