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