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