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