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