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