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