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