Actual source code: options.c
petsc-master 2020-08-25
1: /* Define Feature test macros to make sure atoll is available (SVr4, POSIX.1-2001, 4.3BSD, C99), not in (C89 and POSIX.1-1996) */
2: #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for atoll() */
4: /*
5: These routines simplify the use of command line, file options, etc., and are used to manipulate the options database.
6: This provides the low-level interface, the high level interface is in aoptions.c
8: Some routines use regular malloc and free because it cannot know what malloc is requested with the
9: options database until it has already processed the input.
10: */
12: #include <petsc/private/petscimpl.h>
13: #include <petscviewer.h>
14: #include <ctype.h>
15: #if defined(PETSC_HAVE_MALLOC_H)
16: #include <malloc.h>
17: #endif
18: #if defined(PETSC_HAVE_STRINGS_H)
19: # include <strings.h> /* strcasecmp */
20: #endif
21: #if defined(PETSC_HAVE_YAML)
22: #include <yaml.h>
23: #endif
25: #if defined(PETSC_HAVE_STRCASECMP)
26: #define PetscOptNameCmp(a,b) strcasecmp(a,b)
27: #elif defined(PETSC_HAVE_STRICMP)
28: #define PetscOptNameCmp(a,b) stricmp(a,b)
29: #else
30: #define PetscOptNameCmp(a,b) Error_strcasecmp_not_found
31: #endif
33: #include <petsc/private/hashtable.h>
35: /* This assumes ASCII encoding and ignores locale settings */
36: /* Using tolower() is about 2X slower in microbenchmarks */
37: PETSC_STATIC_INLINE int PetscToLower(int c)
38: {
39: return ((c >= 'A') & (c <= 'Z')) ? c + 'a' - 'A' : c;
40: }
42: /* Bob Jenkins's one at a time hash function (case-insensitive) */
43: PETSC_STATIC_INLINE unsigned int PetscOptHash(const char key[])
44: {
45: unsigned int hash = 0;
46: while (*key) {
47: hash += PetscToLower(*key++);
48: hash += hash << 10;
49: hash ^= hash >> 6;
50: }
51: hash += hash << 3;
52: hash ^= hash >> 11;
53: hash += hash << 15;
54: return hash;
55: }
57: PETSC_STATIC_INLINE int PetscOptEqual(const char a[],const char b[])
58: {
59: return !PetscOptNameCmp(a,b);
60: }
62: KHASH_INIT(HO, kh_cstr_t, int, 1, PetscOptHash, PetscOptEqual)
64: /*
65: This table holds all the options set by the user. For simplicity, we use a static size database
66: */
67: #define MAXOPTNAME 512
68: #define MAXOPTIONS 512
69: #define MAXALIASES 25
70: #define MAXPREFIXES 25
71: #define MAXOPTIONSMONITORS 5
73: struct _n_PetscOptions {
74: PetscOptions previous;
75: int N; /* number of options */
76: char *names[MAXOPTIONS]; /* option names */
77: char *values[MAXOPTIONS]; /* option values */
78: PetscBool used[MAXOPTIONS]; /* flag option use */
79: PetscBool precedentProcessed;
81: /* Hash table */
82: khash_t(HO) *ht;
84: /* Prefixes */
85: int prefixind;
86: int prefixstack[MAXPREFIXES];
87: char prefix[MAXOPTNAME];
89: /* Aliases */
90: int Naliases; /* number or aliases */
91: char *aliases1[MAXALIASES]; /* aliased */
92: char *aliases2[MAXALIASES]; /* aliasee */
94: /* Help */
95: PetscBool help; /* flag whether "-help" is in the database */
96: PetscBool help_intro; /* flag whether "-help intro" is in the database */
98: /* Monitors */
99: PetscBool monitorFromOptions, monitorCancel;
100: PetscErrorCode (*monitor[MAXOPTIONSMONITORS])(const char[],const char[],void*); /* returns control to user after */
101: PetscErrorCode (*monitordestroy[MAXOPTIONSMONITORS])(void**); /* */
102: void *monitorcontext[MAXOPTIONSMONITORS]; /* to pass arbitrary user data into monitor */
103: PetscInt numbermonitors; /* to, for instance, detect options being set */
104: };
106: static PetscOptions defaultoptions = NULL; /* the options database routines query this object for options */
108: /* list of options which preceed others, i.e., are processed in PetscOptionsProcessPrecedentFlags() */
109: static const char *precedentOptions[] = {"-options_monitor","-options_monitor_cancel","-help","-skip_petscrc","-options_file_yaml","-options_string_yaml"};
110: enum PetscPrecedentOption {PO_OPTIONS_MONITOR,PO_OPTIONS_MONITOR_CANCEL,PO_HELP,PO_SKIP_PETSCRC,PO_OPTIONS_FILE_YAML,PO_OPTIONS_STRING_YAML,PO_NUM};
112: static PetscErrorCode PetscOptionsSetValue_Private(PetscOptions,const char[],const char[],int*);
114: /*
115: Options events monitor
116: */
117: static PetscErrorCode PetscOptionsMonitor(PetscOptions options,const char name[],const char value[])
118: {
119: PetscInt i;
122: if (!PetscErrorHandlingInitialized) return 0;
124: if (!value) value = "";
125: if (options->monitorFromOptions) {
126: PetscOptionsMonitorDefault(name,value,NULL);
127: }
128: for (i=0; i<options->numbermonitors; i++) {
129: (*options->monitor[i])(name,value,options->monitorcontext[i]);
130: }
131: return(0);
132: }
134: /*@
135: PetscOptionsCreate - Creates an empty options database.
137: Logically collective
139: Output Parameter:
140: . options - Options database object
142: Level: advanced
144: Developer Note: We may want eventually to pass a MPI_Comm to determine the ownership of the object
146: .seealso: PetscOptionsDestroy(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsInsert(), PetscOptionsSetValue()
147: @*/
148: PetscErrorCode PetscOptionsCreate(PetscOptions *options)
149: {
150: if (!options) return PETSC_ERR_ARG_NULL;
151: *options = (PetscOptions)calloc(1,sizeof(**options));
152: if (!*options) return PETSC_ERR_MEM;
153: return 0;
154: }
156: /*@
157: PetscOptionsDestroy - Destroys an option database.
159: Logically collective on whatever communicator was associated with the call to PetscOptionsCreate()
161: Input Parameter:
162: . options - the PetscOptions object
164: Level: advanced
166: .seealso: PetscOptionsInsert(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsInsert(), PetscOptionsSetValue()
167: @*/
168: PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
169: {
172: if (!*options) return 0;
173: if ((*options)->previous) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You are destroying an option that has been used with PetscOptionsPush() but does not have a corresponding PetscOptionsPop()");
174: PetscOptionsClear(*options);if (ierr) return ierr;
175: /* XXX what about monitors ? */
176: free(*options);
177: *options = NULL;
178: return(0);
179: }
181: /*
182: PetscOptionsCreateDefault - Creates the default global options database
183: */
184: PetscErrorCode PetscOptionsCreateDefault(void)
185: {
188: if (!defaultoptions) {
189: PetscOptionsCreate(&defaultoptions);if (ierr) return ierr;
190: }
191: return 0;
192: }
194: /*@
195: PetscOptionsPush - Push a new PetscOptions object as the default provider of options
196: Allows using different parts of a code to use different options databases
198: Logically Collective
200: Input Parameter:
201: . opt - the options obtained with PetscOptionsCreate()
203: Notes:
204: Use PetscOptionsPop() to return to the previous default options database
206: The collectivity of this routine is complex; only the MPI processes that call this routine will
207: have the affect of these options. If some processes that create objects call this routine and others do
208: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
209: on different ranks.
211: Level: advanced
213: .seealso: PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsInsert(), PetscOptionsSetValue(), PetscOptionsLeft()
215: @*/
216: PetscErrorCode PetscOptionsPush(PetscOptions opt)
217: {
221: PetscOptionsCreateDefault();
222: opt->previous = defaultoptions;
223: defaultoptions = opt;
224: return(0);
225: }
227: /*@
228: PetscOptionsPop - Pop the most recent PetscOptionsPush() to return to the previous default options
230: Logically collective on whatever communicator was associated with the call to PetscOptionsCreate()
232: Notes:
233: Use PetscOptionsPop() to return to the previous default options database
234: Allows using different parts of a code to use different options databases
236: Level: advanced
238: .seealso: PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsInsert(), PetscOptionsSetValue(), PetscOptionsLeft()
240: @*/
241: PetscErrorCode PetscOptionsPop(void)
242: {
243: PetscOptions current = defaultoptions;
246: if (!defaultoptions) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing default options");
247: if (!defaultoptions->previous) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PetscOptionsPop() called too many times");
248: defaultoptions = defaultoptions->previous;
249: current->previous = NULL;
250: return(0);
251: }
253: /*
254: PetscOptionsDestroyDefault - Destroys the default global options database
255: */
256: PetscErrorCode PetscOptionsDestroyDefault(void)
257: {
259: PetscOptions tmp;
261: /* Destroy any options that the user forgot to pop */
262: while (defaultoptions->previous) {
263: tmp = defaultoptions;
264: PetscOptionsPop();
265: PetscOptionsDestroy(&tmp);
266: }
267: PetscOptionsDestroy(&defaultoptions);if (ierr) return ierr;
268: return 0;
269: }
271: /*@C
272: PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.
274: Not collective
276: Input Parameter:
277: . key - string to check if valid
279: Output Parameter:
280: . valid - PETSC_TRUE if a valid key
282: Level: intermediate
283: @*/
284: PetscErrorCode PetscOptionsValidKey(const char key[],PetscBool *valid)
285: {
286: char *ptr;
291: *valid = PETSC_FALSE;
292: if (!key) return(0);
293: if (key[0] != '-') return(0);
294: if (key[1] == '-') key++;
295: if (!isalpha((int)key[1])) return(0);
296: (void) strtod(key,&ptr);
297: if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) return(0);
298: *valid = PETSC_TRUE;
299: return(0);
300: }
302: /*@C
303: PetscOptionsInsertString - Inserts options into the database from a string
305: Logically Collective
307: Input Parameter:
308: + options - options object
309: - in_str - string that contains options separated by blanks
311: Level: intermediate
313: The collectivity of this routine is complex; only the MPI processes that call this routine will
314: have the affect of these options. If some processes that create objects call this routine and others do
315: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
316: on different ranks.
318: Contributed by Boyana Norris
320: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
321: PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
322: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
323: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
324: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
325: PetscOptionsFList(), PetscOptionsEList(), PetscOptionsInsertFile()
326: @*/
327: PetscErrorCode PetscOptionsInsertString(PetscOptions options,const char in_str[])
328: {
329: char *first,*second;
331: PetscToken token;
332: PetscBool key,ispush,ispop,isopts;
335: PetscTokenCreate(in_str,' ',&token);
336: PetscTokenFind(token,&first);
337: while (first) {
338: PetscStrcasecmp(first,"-prefix_push",&ispush);
339: PetscStrcasecmp(first,"-prefix_pop",&ispop);
340: PetscStrcasecmp(first,"-options_file",&isopts);
341: PetscOptionsValidKey(first,&key);
342: if (ispush) {
343: PetscTokenFind(token,&second);
344: PetscOptionsPrefixPush(options,second);
345: PetscTokenFind(token,&first);
346: } else if (ispop) {
347: PetscOptionsPrefixPop(options);
348: PetscTokenFind(token,&first);
349: } else if (isopts) {
350: PetscTokenFind(token,&second);
351: PetscOptionsInsertFile(PETSC_COMM_SELF,options,second,PETSC_TRUE);
352: PetscTokenFind(token,&first);
353: } else if (key) {
354: PetscTokenFind(token,&second);
355: PetscOptionsValidKey(second,&key);
356: if (!key) {
357: PetscOptionsSetValue(options,first,second);
358: PetscTokenFind(token,&first);
359: } else {
360: PetscOptionsSetValue(options,first,NULL);
361: first = second;
362: }
363: } else {
364: PetscTokenFind(token,&first);
365: }
366: }
367: PetscTokenDestroy(&token);
368: return(0);
369: }
371: /*
372: Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
373: */
374: static char *Petscgetline(FILE * f)
375: {
376: size_t size = 0;
377: size_t len = 0;
378: size_t last = 0;
379: char *buf = NULL;
381: if (feof(f)) return NULL;
382: do {
383: size += 1024; /* BUFSIZ is defined as "the optimal read size for this platform" */
384: buf = (char*)realloc((void*)buf,size); /* realloc(NULL,n) is the same as malloc(n) */
385: /* Actually do the read. Note that fgets puts a terminal '\0' on the
386: end of the string, so we make sure we overwrite this */
387: if (!fgets(buf+len,1024,f)) buf[len]=0;
388: PetscStrlen(buf,&len);
389: last = len - 1;
390: } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
391: if (len) return buf;
392: free(buf);
393: return NULL;
394: }
396: /*@C
397: PetscOptionsInsertFile - Inserts options into the database from a file.
399: Collective
401: Input Parameter:
402: + comm - the processes that will share the options (usually PETSC_COMM_WORLD)
403: . options - options database, use NULL for default global database
404: . file - name of file
405: - require - if PETSC_TRUE will generate an error if the file does not exist
408: Notes:
409: Use # for lines that are comments and which should be ignored.
410: Usually, instead of using this command, one should list the file name in the call to PetscInitialize(), this insures that certain options
411: such as -log_view or -malloc_debug are processed properly. This routine only sets options into the options database that will be processed by later
412: calls to XXXSetFromOptions() it should not be used for options listed under PetscInitialize().
413: The collectivity of this routine is complex; only the MPI processes in comm will
414: have the affect of these options. If some processes that create objects call this routine and others do
415: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
416: on different ranks.
418: Level: developer
420: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
421: PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
422: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
423: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
424: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
425: PetscOptionsFList(), PetscOptionsEList()
427: @*/
428: PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm,PetscOptions options,const char file[],PetscBool require)
429: {
430: char *string,fname[PETSC_MAX_PATH_LEN],*vstring = NULL,*astring = NULL,*packed = NULL;
431: char *tokens[4];
433: size_t i,len,bytes;
434: FILE *fd;
435: PetscToken token=NULL;
436: int err;
437: char *cmatch;
438: const char cmt='#';
439: PetscInt line=1;
440: PetscMPIInt rank,cnt=0,acnt=0,counts[2];
441: PetscBool isdir,alias=PETSC_FALSE,valid;
444: MPI_Comm_rank(comm,&rank);
445: PetscMemzero(tokens,sizeof(tokens));
446: if (!rank) {
447: cnt = 0;
448: acnt = 0;
450: PetscFixFilename(file,fname);
451: fd = fopen(fname,"r");
452: PetscTestDirectory(fname,'r',&isdir);
453: if (isdir && require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified options file %s is a directory",fname);
454: if (fd && !isdir) {
455: PetscSegBuffer vseg,aseg;
456: PetscSegBufferCreate(1,4000,&vseg);
457: PetscSegBufferCreate(1,2000,&aseg);
459: /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
460: PetscInfo1(NULL,"Opened options file %s\n",file);
462: while ((string = Petscgetline(fd))) {
463: /* eliminate comments from each line */
464: PetscStrchr(string,cmt,&cmatch);
465: if (cmatch) *cmatch = 0;
466: PetscStrlen(string,&len);
467: /* replace tabs, ^M, \n with " " */
468: for (i=0; i<len; i++) {
469: if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') {
470: string[i] = ' ';
471: }
472: }
473: PetscTokenCreate(string,' ',&token);
474: PetscTokenFind(token,&tokens[0]);
475: if (!tokens[0]) {
476: goto destroy;
477: } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
478: PetscTokenFind(token,&tokens[0]);
479: }
480: for (i=1; i<4; i++) {
481: PetscTokenFind(token,&tokens[i]);
482: }
483: if (!tokens[0]) {
484: goto destroy;
485: } else if (tokens[0][0] == '-') {
486: PetscOptionsValidKey(tokens[0],&valid);
487: if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: invalid option %s",fname,line,tokens[0]);
488: PetscStrlen(tokens[0],&len);
489: PetscSegBufferGet(vseg,len+1,&vstring);
490: PetscArraycpy(vstring,tokens[0],len);
491: vstring[len] = ' ';
492: if (tokens[1]) {
493: PetscOptionsValidKey(tokens[1],&valid);
494: if (valid) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: cannot specify two options per line (%s %s)",fname,line,tokens[0],tokens[1]);
495: PetscStrlen(tokens[1],&len);
496: PetscSegBufferGet(vseg,len+3,&vstring);
497: vstring[0] = '"';
498: PetscArraycpy(vstring+1,tokens[1],len);
499: vstring[len+1] = '"';
500: vstring[len+2] = ' ';
501: }
502: } else {
503: PetscStrcasecmp(tokens[0],"alias",&alias);
504: if (alias) {
505: PetscOptionsValidKey(tokens[1],&valid);
506: if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: invalid aliased option %s",fname,line,tokens[1]);
507: if (!tokens[2]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: alias missing for %s",fname,line,tokens[1]);
508: PetscOptionsValidKey(tokens[2],&valid);
509: if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: invalid aliasee option %s",fname,line,tokens[2]);
510: PetscStrlen(tokens[1],&len);
511: PetscSegBufferGet(aseg,len+1,&astring);
512: PetscArraycpy(astring,tokens[1],len);
513: astring[len] = ' ';
515: PetscStrlen(tokens[2],&len);
516: PetscSegBufferGet(aseg,len+1,&astring);
517: PetscArraycpy(astring,tokens[2],len);
518: astring[len] = ' ';
519: } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown first token in options file %s line %D: %s",fname,line,tokens[0]);
520: }
521: {
522: const char *extraToken = alias ? tokens[3] : tokens[2];
523: if (extraToken) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: extra token %s",fname,line,extraToken);
524: }
525: destroy:
526: free(string);
527: PetscTokenDestroy(&token);
528: alias = PETSC_FALSE;
529: line++;
530: }
531: err = fclose(fd);
532: if (err) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file %s",fname);
533: PetscSegBufferGetSize(aseg,&bytes); /* size without null termination */
534: PetscMPIIntCast(bytes,&acnt);
535: PetscSegBufferGet(aseg,1,&astring);
536: astring[0] = 0;
537: PetscSegBufferGetSize(vseg,&bytes); /* size without null termination */
538: PetscMPIIntCast(bytes,&cnt);
539: PetscSegBufferGet(vseg,1,&vstring);
540: vstring[0] = 0;
541: PetscMalloc1(2+acnt+cnt,&packed);
542: PetscSegBufferExtractTo(aseg,packed);
543: PetscSegBufferExtractTo(vseg,packed+acnt+1);
544: PetscSegBufferDestroy(&aseg);
545: PetscSegBufferDestroy(&vseg);
546: } else if (require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Unable to open options file %s",fname);
547: }
549: counts[0] = acnt;
550: counts[1] = cnt;
551: MPI_Bcast(counts,2,MPI_INT,0,comm);
552: acnt = counts[0];
553: cnt = counts[1];
554: if (rank) {
555: PetscMalloc1(2+acnt+cnt,&packed);
556: }
557: if (acnt || cnt) {
558: MPI_Bcast(packed,2+acnt+cnt,MPI_CHAR,0,comm);
559: astring = packed;
560: vstring = packed + acnt + 1;
561: }
563: if (acnt) {
564: PetscTokenCreate(astring,' ',&token);
565: PetscTokenFind(token,&tokens[0]);
566: while (tokens[0]) {
567: PetscTokenFind(token,&tokens[1]);
568: PetscOptionsSetAlias(options,tokens[0],tokens[1]);
569: PetscTokenFind(token,&tokens[0]);
570: }
571: PetscTokenDestroy(&token);
572: }
574: if (cnt) {
575: PetscOptionsInsertString(options,vstring);
576: }
577: PetscFree(packed);
578: return(0);
579: }
581: static PetscErrorCode PetscOptionsInsertArgs(PetscOptions options,int argc,char *args[])
582: {
584: int left = argc - 1;
585: char **eargs = args + 1;
588: while (left) {
589: PetscBool isoptions_file,isprefixpush,isprefixpop,isp4,tisp4,isp4yourname,isp4rmrank,key;
590: PetscStrcasecmp(eargs[0],"-options_file",&isoptions_file);
591: PetscStrcasecmp(eargs[0],"-prefix_push",&isprefixpush);
592: PetscStrcasecmp(eargs[0],"-prefix_pop",&isprefixpop);
593: PetscStrcasecmp(eargs[0],"-p4pg",&isp4);
594: PetscStrcasecmp(eargs[0],"-p4yourname",&isp4yourname);
595: PetscStrcasecmp(eargs[0],"-p4rmrank",&isp4rmrank);
596: PetscStrcasecmp(eargs[0],"-p4wd",&tisp4);
597: isp4 = (PetscBool) (isp4 || tisp4);
598: PetscStrcasecmp(eargs[0],"-np",&tisp4);
599: isp4 = (PetscBool) (isp4 || tisp4);
600: PetscStrcasecmp(eargs[0],"-p4amslave",&tisp4);
601: PetscOptionsValidKey(eargs[0],&key);
603: if (!key) {
604: eargs++; left--;
605: } else if (isoptions_file) {
606: if (left <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file filename option");
607: if (eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file filename option");
608: PetscOptionsInsertFile(PETSC_COMM_WORLD,options,eargs[1],PETSC_TRUE);
609: eargs += 2; left -= 2;
610: } else if (isprefixpush) {
611: if (left <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option");
612: if (eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option (prefixes cannot start with '-')");
613: PetscOptionsPrefixPush(options,eargs[1]);
614: eargs += 2; left -= 2;
615: } else if (isprefixpop) {
616: PetscOptionsPrefixPop(options);
617: eargs++; left--;
619: /*
620: These are "bad" options that MPICH, etc put on the command line
621: we strip them out here.
622: */
623: } else if (tisp4 || isp4rmrank) {
624: eargs += 1; left -= 1;
625: } else if (isp4 || isp4yourname) {
626: eargs += 2; left -= 2;
627: } else {
628: PetscBool nextiskey = PETSC_FALSE;
629: if (left >= 2) {PetscOptionsValidKey(eargs[1],&nextiskey);}
630: if (left < 2 || nextiskey) {
631: PetscOptionsSetValue(options,eargs[0],NULL);
632: eargs++; left--;
633: } else {
634: PetscOptionsSetValue(options,eargs[0],eargs[1]);
635: eargs += 2; left -= 2;
636: }
637: }
638: }
639: return(0);
640: }
642: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt,const char *val[],PetscBool set[],PetscBool *flg)
643: {
647: if (set[opt]) {
648: PetscOptionsStringToBool(val[opt],flg);
649: } else *flg = PETSC_FALSE;
650: return(0);
651: }
653: /* Process options with absolute precedence */
654: static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options,int argc,char *args[],PetscBool *skip_petscrc,PetscBool *skip_petscrc_set)
655: {
656: const char* const *opt = precedentOptions;
657: const size_t n = PO_NUM;
658: size_t o;
659: int a;
660: const char **val;
661: PetscBool *set;
662: PetscErrorCode ierr;
665: PetscCalloc2(n,&val,n,&set);
667: /* Look for options possibly set using PetscOptionsSetValue beforehand */
668: for (o=0; o<n; o++) {
669: PetscOptionsFindPair(options,NULL,opt[o],&val[o],&set[o]);
670: }
672: /* Loop through all args to collect last occuring value of each option */
673: for (a=1; a<argc; a++) {
674: PetscBool valid, eq;
676: PetscOptionsValidKey(args[a],&valid);
677: if (!valid) continue;
678: for (o=0; o<n; o++) {
679: PetscStrcasecmp(args[a],opt[o],&eq);
680: if (eq) {
681: set[o] = PETSC_TRUE;
682: if (a == argc-1 || !args[a+1] || !args[a+1][0] || args[a+1][0] == '-') val[o] = NULL;
683: else val[o] = args[a+1];
684: break;
685: }
686: }
687: }
689: /* Process flags */
690: PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro);
691: if (options->help_intro) options->help = PETSC_TRUE;
692: else {PetscOptionsStringToBoolIfSet_Private(PO_HELP, val,set,&options->help);}
693: PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL,val,set,&options->monitorCancel);
694: PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR, val,set,&options->monitorFromOptions);
695: PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC, val,set,skip_petscrc);
696: *skip_petscrc_set = set[PO_SKIP_PETSCRC];
698: /* Store precedent options in database and mark them as used */
699: for (o=0; o<n; o++) {
700: if (set[o]) {
701: int pos;
703: PetscOptionsSetValue_Private(options,opt[o],val[o],&pos);
704: options->used[pos] = PETSC_TRUE;
705: }
706: }
708: PetscFree2(val,set);
709: options->precedentProcessed = PETSC_TRUE;
710: return(0);
711: }
713: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options,const char name[],PetscBool *flg)
714: {
715: int i;
718: *flg = PETSC_FALSE;
719: if (options->precedentProcessed) {
720: for (i=0; i<PO_NUM; i++) {
721: if (!PetscOptNameCmp(precedentOptions[i],name)) {
722: /* check if precedent option has been set already */
723: PetscOptionsFindPair(options,NULL,name,NULL,flg);
724: if (*flg) break;
725: }
726: }
727: }
728: return(0);
729: }
731: /*@C
732: PetscOptionsInsert - Inserts into the options database from the command line,
733: the environmental variable and a file.
735: Collective on PETSC_COMM_WORLD
737: Input Parameters:
738: + options - options database or NULL for the default global database
739: . argc - count of number of command line arguments
740: . args - the command line arguments
741: - file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
742: Use NULL to not check for code specific file.
743: Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.
745: Note:
746: Since PetscOptionsInsert() is automatically called by PetscInitialize(),
747: the user does not typically need to call this routine. PetscOptionsInsert()
748: can be called several times, adding additional entries into the database.
750: Options Database Keys:
751: . -options_file <filename> - read options from a file
753: See PetscInitialize() for options related to option database monitoring.
755: Level: advanced
757: .seealso: PetscOptionsDestroy(), PetscOptionsView(), PetscOptionsInsertString(), PetscOptionsInsertFile(),
758: PetscInitialize()
759: @*/
760: PetscErrorCode PetscOptionsInsert(PetscOptions options,int *argc,char ***args,const char file[])
761: {
763: PetscMPIInt rank;
764: char filename[PETSC_MAX_PATH_LEN];
765: PetscBool hasArgs = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
766: PetscBool skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;
770: if (hasArgs && !(args && *args)) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
771: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
773: if (!options) {
774: PetscOptionsCreateDefault();
775: options = defaultoptions;
776: }
777: if (hasArgs) {
778: /* process options with absolute precedence */
779: PetscOptionsProcessPrecedentFlags(options,*argc,*args,&skipPetscrc,&skipPetscrcSet);
780: }
781: if (file && file[0]) {
782: PetscStrreplace(PETSC_COMM_WORLD,file,filename,PETSC_MAX_PATH_LEN);
783: PetscOptionsInsertFile(PETSC_COMM_WORLD,options,filename,PETSC_TRUE);
784: /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
785: if (!skipPetscrcSet) {PetscOptionsGetBool(options,NULL,"-skip_petscrc",&skipPetscrc,NULL);}
786: }
787: if (!skipPetscrc) {
788: PetscGetHomeDirectory(filename,PETSC_MAX_PATH_LEN-16);
789: /* PetscOptionsInsertFile() does a fopen() on rank0 only - so only rank0 HomeDir value is relavent */
790: if (filename[0]) { PetscStrcat(filename,"/.petscrc"); }
791: PetscOptionsInsertFile(PETSC_COMM_WORLD,options,filename,PETSC_FALSE);
792: PetscOptionsInsertFile(PETSC_COMM_WORLD,options,".petscrc",PETSC_FALSE);
793: PetscOptionsInsertFile(PETSC_COMM_WORLD,options,"petscrc",PETSC_FALSE);
794: }
796: /* insert environment options */
797: {
798: char *eoptions = NULL;
799: size_t len = 0;
800: if (!rank) {
801: eoptions = (char*)getenv("PETSC_OPTIONS");
802: PetscStrlen(eoptions,&len);
803: MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
804: } else {
805: MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
806: if (len) {
807: PetscMalloc1(len+1,&eoptions);
808: }
809: }
810: if (len) {
811: MPI_Bcast(eoptions,len,MPI_CHAR,0,PETSC_COMM_WORLD);
812: if (rank) eoptions[len] = 0;
813: PetscOptionsInsertString(options,eoptions);
814: if (rank) {PetscFree(eoptions);}
815: }
816: }
818: #if defined(PETSC_HAVE_YAML)
819: {
820: char *eoptions = NULL;
821: size_t len = 0;
822: if (!rank) {
823: eoptions = (char*)getenv("PETSC_OPTIONS_YAML");
824: PetscStrlen(eoptions,&len);
825: MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
826: } else {
827: MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
828: if (len) {
829: PetscMalloc1(len+1,&eoptions);
830: }
831: }
832: if (len) {
833: MPI_Bcast(eoptions,len,MPI_CHAR,0,PETSC_COMM_WORLD);
834: if (rank) eoptions[len] = 0;
835: PetscOptionsInsertStringYAML(options,eoptions);
836: if (rank) {PetscFree(eoptions);}
837: }
838: }
839: {
840: char yaml_file[PETSC_MAX_PATH_LEN];
841: char yaml_string[BUFSIZ];
842: PetscBool yaml_flg;
843: PetscOptionsGetString(NULL,NULL,"-options_file_yaml",yaml_file,sizeof(yaml_file),&yaml_flg);
844: if (yaml_flg) {
845: PetscOptionsInsertFileYAML(PETSC_COMM_WORLD,yaml_file,PETSC_TRUE);
846: }
847: PetscOptionsGetString(NULL,NULL,"-options_string_yaml",yaml_string,sizeof(yaml_string),&yaml_flg);
848: if (yaml_flg) {
849: PetscOptionsInsertStringYAML(NULL,yaml_string);
850: }
851: }
852: #endif
854: /* insert command line options here because they take precedence over arguments in petscrc/environment */
855: if (hasArgs) {PetscOptionsInsertArgs(options,*argc,*args);}
856: return(0);
857: }
859: /*@C
860: PetscOptionsView - Prints the options that have been loaded. This is
861: useful for debugging purposes.
863: Logically Collective on PetscViewer
865: Input Parameter:
866: + options - options database, use NULL for default global database
867: - viewer - must be an PETSCVIEWERASCII viewer
869: Options Database Key:
870: . -options_view - Activates PetscOptionsView() within PetscFinalize()
872: Notes:
873: Only the rank zero process of MPI_Comm used to create view prints the option values. Other processes
874: may have different values but they are not printed.
876: Level: advanced
878: .seealso: PetscOptionsAllUsed()
879: @*/
880: PetscErrorCode PetscOptionsView(PetscOptions options,PetscViewer viewer)
881: {
883: PetscInt i;
884: PetscBool isascii;
888: options = options ? options : defaultoptions;
889: if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
890: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
891: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Only supports ASCII viewer");
893: if (!options->N) {
894: PetscViewerASCIIPrintf(viewer,"#No PETSc Option Table entries\n");
895: return(0);
896: }
898: PetscViewerASCIIPrintf(viewer,"#PETSc Option Table entries:\n");
899: for (i=0; i<options->N; i++) {
900: if (options->values[i]) {
901: PetscViewerASCIIPrintf(viewer,"-%s %s\n",options->names[i],options->values[i]);
902: } else {
903: PetscViewerASCIIPrintf(viewer,"-%s\n",options->names[i]);
904: }
905: }
906: PetscViewerASCIIPrintf(viewer,"#End of PETSc Option Table entries\n");
907: return(0);
908: }
910: /*
911: Called by error handlers to print options used in run
912: */
913: PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
914: {
915: PetscInt i;
916: PetscOptions options = defaultoptions;
919: if (options->N) {
920: (*PetscErrorPrintf)("PETSc Option Table entries:\n");
921: } else {
922: (*PetscErrorPrintf)("No PETSc Option Table entries\n");
923: }
924: for (i=0; i<options->N; i++) {
925: if (options->values[i]) {
926: (*PetscErrorPrintf)("-%s %s\n",options->names[i],options->values[i]);
927: } else {
928: (*PetscErrorPrintf)("-%s\n",options->names[i]);
929: }
930: }
931: return(0);
932: }
934: /*@C
935: PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.
937: Logically Collective
939: Input Parameter:
940: + options - options database, or NULL for the default global database
941: - prefix - The string to append to the existing prefix
943: Options Database Keys:
944: + -prefix_push <some_prefix_> - push the given prefix
945: - -prefix_pop - pop the last prefix
947: Notes:
948: It is common to use this in conjunction with -options_file as in
950: $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop
952: where the files no longer require all options to be prefixed with -system2_.
954: The collectivity of this routine is complex; only the MPI processes that call this routine will
955: have the affect of these options. If some processes that create objects call this routine and others do
956: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
957: on different ranks.
959: Level: advanced
961: .seealso: PetscOptionsPrefixPop(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
962: @*/
963: PetscErrorCode PetscOptionsPrefixPush(PetscOptions options,const char prefix[])
964: {
966: size_t n;
967: PetscInt start;
968: char key[MAXOPTNAME+1];
969: PetscBool valid;
973: options = options ? options : defaultoptions;
974: 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);
975: key[0] = '-'; /* keys must start with '-' */
976: PetscStrncpy(key+1,prefix,sizeof(key)-1);
977: PetscOptionsValidKey(key,&valid);
978: if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Given prefix \"%s\" not valid (the first character must be a letter, do not include leading '-')",prefix);
979: start = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
980: PetscStrlen(prefix,&n);
981: if (n+1 > sizeof(options->prefix)-start) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Maximum prefix length %d exceeded",sizeof(options->prefix));
982: PetscArraycpy(options->prefix+start,prefix,n+1);
983: options->prefixstack[options->prefixind++] = start+n;
984: return(0);
985: }
987: /*@C
988: PetscOptionsPrefixPop - Remove the latest options prefix, see PetscOptionsPrefixPush() for details
990: Logically Collective on the MPI_Comm that called PetscOptionsPrefixPush()
992: Input Parameters:
993: . options - options database, or NULL for the default global database
995: Level: advanced
997: .seealso: PetscOptionsPrefixPush(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
998: @*/
999: PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1000: {
1001: PetscInt offset;
1004: options = options ? options : defaultoptions;
1005: if (options->prefixind < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"More prefixes popped than pushed");
1006: options->prefixind--;
1007: offset = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
1008: options->prefix[offset] = 0;
1009: return(0);
1010: }
1012: /*@C
1013: PetscOptionsClear - Removes all options form the database leaving it empty.
1015: Logically Collective
1017: Input Parameters:
1018: . options - options database, use NULL for the default global database
1020: The collectivity of this routine is complex; only the MPI processes that call this routine will
1021: have the affect of these options. If some processes that create objects call this routine and others do
1022: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1023: on different ranks.
1025: Level: developer
1027: .seealso: PetscOptionsInsert()
1028: @*/
1029: PetscErrorCode PetscOptionsClear(PetscOptions options)
1030: {
1031: PetscInt i;
1033: options = options ? options : defaultoptions;
1034: if (!options) return 0;
1036: for (i=0; i<options->N; i++) {
1037: if (options->names[i]) free(options->names[i]);
1038: if (options->values[i]) free(options->values[i]);
1039: }
1040: options->N = 0;
1042: for (i=0; i<options->Naliases; i++) {
1043: free(options->aliases1[i]);
1044: free(options->aliases2[i]);
1045: }
1046: options->Naliases = 0;
1048: /* destroy hash table */
1049: kh_destroy(HO,options->ht);
1050: options->ht = NULL;
1052: options->prefixind = 0;
1053: options->prefix[0] = 0;
1054: options->help = PETSC_FALSE;
1055: return 0;
1056: }
1058: /*@C
1059: PetscOptionsSetAlias - Makes a key and alias for another key
1061: Logically Collective
1063: Input Parameters:
1064: + options - options database, or NULL for default global database
1065: . newname - the alias
1066: - oldname - the name that alias will refer to
1068: Level: advanced
1070: The collectivity of this routine is complex; only the MPI processes that call this routine will
1071: have the affect of these options. If some processes that create objects call this routine and others do
1072: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1073: on different ranks.
1075: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1076: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
1077: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1078: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1079: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1080: PetscOptionsFList(), PetscOptionsEList()
1081: @*/
1082: PetscErrorCode PetscOptionsSetAlias(PetscOptions options,const char newname[],const char oldname[])
1083: {
1084: PetscInt n;
1085: size_t len;
1086: PetscBool valid;
1092: options = options ? options : defaultoptions;
1093: PetscOptionsValidKey(newname,&valid);
1094: if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliased option %s",newname);
1095: PetscOptionsValidKey(oldname,&valid);
1096: if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliasee option %s",oldname);
1098: n = options->Naliases;
1099: 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);
1101: newname++; oldname++;
1102: PetscStrlen(newname,&len);
1103: options->aliases1[n] = (char*)malloc((len+1)*sizeof(char));
1104: PetscStrcpy(options->aliases1[n],newname);
1105: PetscStrlen(oldname,&len);
1106: options->aliases2[n] = (char*)malloc((len+1)*sizeof(char));
1107: PetscStrcpy(options->aliases2[n],oldname);
1108: options->Naliases++;
1109: return(0);
1110: }
1112: /*@C
1113: PetscOptionsSetValue - Sets an option name-value pair in the options
1114: database, overriding whatever is already present.
1116: Logically Collective
1118: Input Parameters:
1119: + options - options database, use NULL for the default global database
1120: . name - name of option, this SHOULD have the - prepended
1121: - value - the option value (not used for all options, so can be NULL)
1123: Level: intermediate
1125: Note:
1126: This function can be called BEFORE PetscInitialize()
1128: The collectivity of this routine is complex; only the MPI processes that call this routine will
1129: have the affect of these options. If some processes that create objects call this routine and others do
1130: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1131: on different ranks.
1133: Developers Note: Uses malloc() directly because PETSc may not be initialized yet.
1135: .seealso: PetscOptionsInsert(), PetscOptionsClearValue()
1136: @*/
1137: PetscErrorCode PetscOptionsSetValue(PetscOptions options,const char name[],const char value[])
1138: {
1139: return PetscOptionsSetValue_Private(options,name,value,NULL);
1140: }
1142: static PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options,const char name[],const char value[],int *pos)
1143: {
1144: size_t len;
1145: int N,n,i;
1146: char **names;
1147: char fullname[MAXOPTNAME] = "";
1148: PetscBool flg;
1151: if (!options) {
1152: PetscOptionsCreateDefault();if (ierr) return ierr;
1153: options = defaultoptions;
1154: }
1156: if (name[0] != '-') return PETSC_ERR_ARG_OUTOFRANGE;
1158: PetscOptionsSkipPrecedent(options,name,&flg);
1159: if (flg) return 0;
1161: name++; /* skip starting dash */
1163: if (options->prefixind > 0) {
1164: strncpy(fullname,options->prefix,sizeof(fullname));
1165: fullname[sizeof(fullname)-1] = 0;
1166: strncat(fullname,name,sizeof(fullname)-strlen(fullname)-1);
1167: fullname[sizeof(fullname)-1] = 0;
1168: name = fullname;
1169: }
1171: /* check against aliases */
1172: N = options->Naliases;
1173: for (i=0; i<N; i++) {
1174: int result = PetscOptNameCmp(options->aliases1[i],name);
1175: if (!result) { name = options->aliases2[i]; break; }
1176: }
1178: /* slow search */
1179: N = n = options->N;
1180: names = options->names;
1181: for (i=0; i<N; i++) {
1182: int result = PetscOptNameCmp(names[i],name);
1183: if (!result) {
1184: n = i; goto setvalue;
1185: } else if (result > 0) {
1186: n = i; break;
1187: }
1188: }
1189: if (N >= MAXOPTIONS) return PETSC_ERR_MEM;
1190: /* shift remaining values up 1 */
1191: for (i=N; i>n; i--) {
1192: options->names[i] = options->names[i-1];
1193: options->values[i] = options->values[i-1];
1194: options->used[i] = options->used[i-1];
1195: }
1196: options->names[n] = NULL;
1197: options->values[n] = NULL;
1198: options->used[n] = PETSC_FALSE;
1199: options->N++;
1201: /* destroy hash table */
1202: kh_destroy(HO,options->ht);
1203: options->ht = NULL;
1205: /* set new name */
1206: len = strlen(name);
1207: options->names[n] = (char*)malloc((len+1)*sizeof(char));
1208: if (!options->names[n]) return PETSC_ERR_MEM;
1209: strcpy(options->names[n],name);
1211: setvalue:
1212: /* set new value */
1213: if (options->values[n]) free(options->values[n]);
1214: len = value ? strlen(value) : 0;
1215: if (len) {
1216: options->values[n] = (char*)malloc((len+1)*sizeof(char));
1217: if (!options->values[n]) return PETSC_ERR_MEM;
1218: strcpy(options->values[n],value);
1219: } else {
1220: options->values[n] = NULL;
1221: }
1223: if (PetscErrorHandlingInitialized) {
1224: PetscOptionsMonitor(options,name,value);
1225: }
1226: if (pos) *pos = n;
1227: return 0;
1228: }
1230: /*@C
1231: PetscOptionsClearValue - Clears an option name-value pair in the options
1232: database, overriding whatever is already present.
1234: Logically Collective
1236: Input Parameter:
1237: + options - options database, use NULL for the default global database
1238: - name - name of option, this SHOULD have the - prepended
1240: Level: intermediate
1242: The collectivity of this routine is complex; only the MPI processes that call this routine will
1243: have the affect of these options. If some processes that create objects call this routine and others do
1244: not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1245: on different ranks.
1247: .seealso: PetscOptionsInsert()
1248: @*/
1249: PetscErrorCode PetscOptionsClearValue(PetscOptions options,const char name[])
1250: {
1251: int N,n,i;
1252: char **names;
1256: options = options ? options : defaultoptions;
1257: if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);
1259: if (!PetscOptNameCmp(name,"-help")) options->help = PETSC_FALSE;
1261: name++; /* skip starting dash */
1263: /* slow search */
1264: N = n = options->N;
1265: names = options->names;
1266: for (i=0; i<N; i++) {
1267: int result = PetscOptNameCmp(names[i],name);
1268: if (!result) {
1269: n = i; break;
1270: } else if (result > 0) {
1271: n = N; break;
1272: }
1273: }
1274: if (n == N) return(0); /* it was not present */
1276: /* remove name and value */
1277: if (options->names[n]) free(options->names[n]);
1278: if (options->values[n]) free(options->values[n]);
1279: /* shift remaining values down 1 */
1280: for (i=n; i<N-1; i++) {
1281: options->names[i] = options->names[i+1];
1282: options->values[i] = options->values[i+1];
1283: options->used[i] = options->used[i+1];
1284: }
1285: options->N--;
1287: /* destroy hash table */
1288: kh_destroy(HO,options->ht);
1289: options->ht = NULL;
1291: PetscOptionsMonitor(options,name,NULL);
1292: return(0);
1293: }
1295: /*@C
1296: PetscOptionsFindPair - Gets an option name-value pair from the options database.
1298: Not Collective
1300: Input Parameters:
1301: + options - options database, use NULL for the default global database
1302: . pre - the string to prepend to the name or NULL, this SHOULD NOT have the "-" prepended
1303: - name - name of option, this SHOULD have the "-" prepended
1305: Output Parameters:
1306: + value - the option value (optional, not used for all options)
1307: - set - whether the option is set (optional)
1309: Notes:
1310: Each process may find different values or no value depending on how options were inserted into the database
1312: Level: developer
1314: .seealso: PetscOptionsSetValue(), PetscOptionsClearValue()
1315: @*/
1316: PetscErrorCode PetscOptionsFindPair(PetscOptions options,const char pre[],const char name[],const char *value[],PetscBool *set)
1317: {
1318: char buf[MAXOPTNAME];
1319: PetscBool usehashtable = PETSC_TRUE;
1320: PetscBool matchnumbers = PETSC_TRUE;
1324: options = options ? options : defaultoptions;
1325: if (pre && PetscUnlikely(pre[0] == '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre);
1326: if (PetscUnlikely(name[0] != '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);
1328: name++; /* skip starting dash */
1330: /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1331: if (pre && pre[0]) {
1332: char *ptr = buf;
1333: if (name[0] == '-') { *ptr++ = '-'; name++; }
1334: PetscStrncpy(ptr,pre,buf+sizeof(buf)-ptr);
1335: PetscStrlcat(buf,name,sizeof(buf));
1336: name = buf;
1337: }
1339: if (PetscDefined(USE_DEBUG)) {
1340: PetscBool valid;
1341: char key[MAXOPTNAME+1] = "-";
1342: PetscStrncpy(key+1,name,sizeof(key)-1);
1343: PetscOptionsValidKey(key,&valid);
1344: if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1345: }
1347: if (!options->ht && usehashtable) {
1348: int i,ret;
1349: khiter_t it;
1350: khash_t(HO) *ht;
1351: ht = kh_init(HO);
1352: if (PetscUnlikely(!ht)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1353: ret = kh_resize(HO,ht,options->N*2); /* twice the required size to reduce risk of collisions */
1354: if (PetscUnlikely(ret)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1355: for (i=0; i<options->N; i++) {
1356: it = kh_put(HO,ht,options->names[i],&ret);
1357: if (PetscUnlikely(ret != 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1358: kh_val(ht,it) = i;
1359: }
1360: options->ht = ht;
1361: }
1363: if (usehashtable)
1364: { /* fast search */
1365: khash_t(HO) *ht = options->ht;
1366: khiter_t it = kh_get(HO,ht,name);
1367: if (it != kh_end(ht)) {
1368: int i = kh_val(ht,it);
1369: options->used[i] = PETSC_TRUE;
1370: if (value) *value = options->values[i];
1371: if (set) *set = PETSC_TRUE;
1372: return(0);
1373: }
1374: } else
1375: { /* slow search */
1376: int i, N = options->N;
1377: for (i=0; i<N; i++) {
1378: int result = PetscOptNameCmp(options->names[i],name);
1379: if (!result) {
1380: options->used[i] = PETSC_TRUE;
1381: if (value) *value = options->values[i];
1382: if (set) *set = PETSC_TRUE;
1383: return(0);
1384: } else if (result > 0) {
1385: break;
1386: }
1387: }
1388: }
1390: /*
1391: The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
1392: Maybe this special lookup mode should be enabled on request with a push/pop API.
1393: The feature of matching _%d_ used sparingly in the codebase.
1394: */
1395: if (matchnumbers) {
1396: int i,j,cnt = 0,locs[16],loce[16];
1397: /* determine the location and number of all _%d_ in the key */
1398: for (i=0; name[i]; i++) {
1399: if (name[i] == '_') {
1400: for (j=i+1; name[j]; j++) {
1401: if (name[j] >= '0' && name[j] <= '9') continue;
1402: if (name[j] == '_' && j > i+1) { /* found a number */
1403: locs[cnt] = i+1;
1404: loce[cnt++] = j+1;
1405: }
1406: i = j-1;
1407: break;
1408: }
1409: }
1410: }
1411: for (i=0; i<cnt; i++) {
1412: PetscBool found;
1413: char opt[MAXOPTNAME+1] = "-", tmp[MAXOPTNAME];
1414: PetscStrncpy(tmp,name,PetscMin((size_t)(locs[i]+1),sizeof(tmp)));
1415: PetscStrlcat(opt,tmp,sizeof(opt));
1416: PetscStrlcat(opt,name+loce[i],sizeof(opt));
1417: PetscOptionsFindPair(options,NULL,opt,value,&found);
1418: if (found) {if (set) *set = PETSC_TRUE; return(0);}
1419: }
1420: }
1422: if (set) *set = PETSC_FALSE;
1423: return(0);
1424: }
1426: /* Check whether any option begins with pre+name */
1427: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[],const char *value[],PetscBool *set)
1428: {
1429: char buf[MAXOPTNAME];
1430: int numCnt = 0, locs[16],loce[16];
1434: options = options ? options : defaultoptions;
1435: if (pre && pre[0] == '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre);
1436: if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);
1438: name++; /* skip starting dash */
1440: /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1441: if (pre && pre[0]) {
1442: char *ptr = buf;
1443: if (name[0] == '-') { *ptr++ = '-'; name++; }
1444: PetscStrncpy(ptr,pre,sizeof(buf)+(size_t)(ptr-buf));
1445: PetscStrlcat(buf,name,sizeof(buf));
1446: name = buf;
1447: }
1449: if (PetscDefined(USE_DEBUG)) {
1450: PetscBool valid;
1451: char key[MAXOPTNAME+1] = "-";
1452: PetscStrncpy(key+1,name,sizeof(key)-1);
1453: PetscOptionsValidKey(key,&valid);
1454: if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1455: }
1457: /* determine the location and number of all _%d_ in the key */
1458: {
1459: int i,j;
1460: for (i=0; name[i]; i++) {
1461: if (name[i] == '_') {
1462: for (j=i+1; name[j]; j++) {
1463: if (name[j] >= '0' && name[j] <= '9') continue;
1464: if (name[j] == '_' && j > i+1) { /* found a number */
1465: locs[numCnt] = i+1;
1466: loce[numCnt++] = j+1;
1467: }
1468: i = j-1;
1469: break;
1470: }
1471: }
1472: }
1473: }
1475: { /* slow search */
1476: int c, i;
1477: size_t len;
1478: PetscBool match;
1480: for (c = -1; c < numCnt; ++c) {
1481: char opt[MAXOPTNAME+1] = "", tmp[MAXOPTNAME];
1483: if (c < 0) {
1484: PetscStrcpy(opt,name);
1485: } else {
1486: PetscStrncpy(tmp,name,PetscMin((size_t)(locs[c]+1),sizeof(tmp)));
1487: PetscStrlcat(opt,tmp,sizeof(opt));
1488: PetscStrlcat(opt,name+loce[c],sizeof(opt));
1489: }
1490: PetscStrlen(opt,&len);
1491: for (i=0; i<options->N; i++) {
1492: PetscStrncmp(options->names[i],opt,len,&match);
1493: if (match) {
1494: options->used[i] = PETSC_TRUE;
1495: if (value) *value = options->values[i];
1496: if (set) *set = PETSC_TRUE;
1497: return(0);
1498: }
1499: }
1500: }
1501: }
1503: if (set) *set = PETSC_FALSE;
1504: return(0);
1505: }
1507: /*@C
1508: PetscOptionsReject - Generates an error if a certain option is given.
1510: Not Collective
1512: Input Parameters:
1513: + options - options database, use NULL for default global database
1514: . pre - the option prefix (may be NULL)
1515: . name - the option name one is seeking
1516: - mess - error message (may be NULL)
1518: Level: advanced
1520: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1521: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1522: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1523: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1524: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1525: PetscOptionsFList(), PetscOptionsEList()
1526: @*/
1527: PetscErrorCode PetscOptionsReject(PetscOptions options,const char pre[],const char name[],const char mess[])
1528: {
1530: PetscBool flag = PETSC_FALSE;
1533: PetscOptionsHasName(options,pre,name,&flag);
1534: if (flag) {
1535: if (mess && mess[0]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s with %s",pre?pre:"",name+1,mess);
1536: else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s",pre?pre:"",name+1);
1537: }
1538: return(0);
1539: }
1541: /*@C
1542: PetscOptionsHasHelp - Determines whether the "-help" option is in the database.
1544: Not Collective
1546: Input Parameters:
1547: . options - options database, use NULL for default global database
1549: Output Parameters:
1550: . set - PETSC_TRUE if found else PETSC_FALSE.
1552: Level: advanced
1554: .seealso: PetscOptionsHasName()
1555: @*/
1556: PetscErrorCode PetscOptionsHasHelp(PetscOptions options,PetscBool *set)
1557: {
1560: options = options ? options : defaultoptions;
1561: *set = options->help;
1562: return(0);
1563: }
1565: PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options,PetscBool *set)
1566: {
1569: options = options ? options : defaultoptions;
1570: *set = options->help_intro;
1571: return(0);
1572: }
1574: /*@C
1575: PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or boolean, even
1576: its value is set to false.
1578: Not Collective
1580: Input Parameters:
1581: + options - options database, use NULL for default global database
1582: . pre - string to prepend to the name or NULL
1583: - name - the option one is seeking
1585: Output Parameters:
1586: . set - PETSC_TRUE if found else PETSC_FALSE.
1588: Level: beginner
1590: Notes:
1591: In many cases you probably want to use PetscOptionsGetBool() instead of calling this, to allowing toggling values.
1593: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1594: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1595: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1596: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1597: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1598: PetscOptionsFList(), PetscOptionsEList()
1599: @*/
1600: PetscErrorCode PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool *set)
1601: {
1602: const char *value;
1604: PetscBool flag;
1607: PetscOptionsFindPair(options,pre,name,&value,&flag);
1608: if (set) *set = flag;
1609: return(0);
1610: }
1612: /*@C
1613: PetscOptionsGetAll - Lists all the options the program was run with in a single string.
1615: Not Collective
1617: Input Parameter:
1618: . options - the options database, use NULL for the default global database
1620: Output Parameter:
1621: . copts - pointer where string pointer is stored
1623: Notes:
1624: The array and each entry in the array should be freed with PetscFree()
1625: Each process may have different values depending on how the options were inserted into the database
1627: Level: advanced
1629: .seealso: PetscOptionsAllUsed(), PetscOptionsView(), PetscOptionsPush(), PetscOptionsPop()
1630: @*/
1631: PetscErrorCode PetscOptionsGetAll(PetscOptions options,char *copts[])
1632: {
1634: PetscInt i;
1635: size_t len = 1,lent = 0;
1636: char *coptions = NULL;
1640: options = options ? options : defaultoptions;
1641: /* count the length of the required string */
1642: for (i=0; i<options->N; i++) {
1643: PetscStrlen(options->names[i],&lent);
1644: len += 2 + lent;
1645: if (options->values[i]) {
1646: PetscStrlen(options->values[i],&lent);
1647: len += 1 + lent;
1648: }
1649: }
1650: PetscMalloc1(len,&coptions);
1651: coptions[0] = 0;
1652: for (i=0; i<options->N; i++) {
1653: PetscStrcat(coptions,"-");
1654: PetscStrcat(coptions,options->names[i]);
1655: PetscStrcat(coptions," ");
1656: if (options->values[i]) {
1657: PetscStrcat(coptions,options->values[i]);
1658: PetscStrcat(coptions," ");
1659: }
1660: }
1661: *copts = coptions;
1662: return(0);
1663: }
1665: /*@C
1666: PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database
1668: Not Collective
1670: Input Parameter:
1671: + options - options database, use NULL for default global database
1672: - name - string name of option
1674: Output Parameter:
1675: . used - PETSC_TRUE if the option was used, otherwise false, including if option was not found in options database
1677: Level: advanced
1679: Notes:
1680: The value returned may be different on each process and depends on which options have been processed
1681: on the given process
1683: .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed()
1684: @*/
1685: PetscErrorCode PetscOptionsUsed(PetscOptions options,const char *name,PetscBool *used)
1686: {
1687: PetscInt i;
1693: options = options ? options : defaultoptions;
1694: *used = PETSC_FALSE;
1695: for (i=0; i<options->N; i++) {
1696: PetscStrcmp(options->names[i],name,used);
1697: if (*used) {
1698: *used = options->used[i];
1699: break;
1700: }
1701: }
1702: return(0);
1703: }
1705: /*@
1706: PetscOptionsAllUsed - Returns a count of the number of options in the
1707: database that have never been selected.
1709: Not Collective
1711: Input Parameter:
1712: . options - options database, use NULL for default global database
1714: Output Parameter:
1715: . N - count of options not used
1717: Level: advanced
1719: Notes:
1720: The value returned may be different on each process and depends on which options have been processed
1721: on the given process
1723: .seealso: PetscOptionsView()
1724: @*/
1725: PetscErrorCode PetscOptionsAllUsed(PetscOptions options,PetscInt *N)
1726: {
1727: PetscInt i,n = 0;
1731: options = options ? options : defaultoptions;
1732: for (i=0; i<options->N; i++) {
1733: if (!options->used[i]) n++;
1734: }
1735: *N = n;
1736: return(0);
1737: }
1739: /*@
1740: PetscOptionsLeft - Prints to screen any options that were set and never used.
1742: Not Collective
1744: Input Parameter:
1745: . options - options database; use NULL for default global database
1747: Options Database Key:
1748: . -options_left - activates PetscOptionsAllUsed() within PetscFinalize()
1750: Notes:
1751: This is rarely used directly, it is called by PetscFinalize() in debug more or if -options_left
1752: is passed otherwise to help users determine possible mistakes in their usage of options. This
1753: only prints values on process zero of PETSC_COMM_WORLD. Other processes depending the objects
1754: used may have different options that are left unused.
1756: Level: advanced
1758: .seealso: PetscOptionsAllUsed()
1759: @*/
1760: PetscErrorCode PetscOptionsLeft(PetscOptions options)
1761: {
1763: PetscInt i;
1764: PetscInt cnt = 0;
1765: PetscOptions toptions;
1768: toptions = options ? options : defaultoptions;
1769: for (i=0; i<toptions->N; i++) {
1770: if (!toptions->used[i]) {
1771: if (toptions->values[i]) {
1772: PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",toptions->names[i],toptions->values[i]);
1773: } else {
1774: PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",toptions->names[i]);
1775: }
1776: }
1777: }
1778: if (!options) {
1779: toptions = defaultoptions;
1780: while (toptions->previous) {
1781: cnt++;
1782: toptions = toptions->previous;
1783: }
1784: if (cnt) {
1785: PetscPrintf(PETSC_COMM_WORLD,"Option left: You may have forgotten some calls to PetscOptionsPop(),\n PetscOptionsPop() has been called %D less times than PetscOptionsPush()\n",cnt);
1786: }
1787: }
1788: return(0);
1789: }
1791: /*@C
1792: PetscOptionsLeftGet - Returns all options that were set and never used.
1794: Not Collective
1796: Input Parameter:
1797: . options - options database, use NULL for default global database
1799: Output Parameter:
1800: + N - count of options not used
1801: . names - names of options not used
1802: - values - values of options not used
1804: Level: advanced
1806: Notes:
1807: Users should call PetscOptionsLeftRestore() to free the memory allocated in this routine
1808: Notes: The value returned may be different on each process and depends on which options have been processed
1809: on the given process
1811: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft()
1812: @*/
1813: PetscErrorCode PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[])
1814: {
1816: PetscInt i,n;
1822: options = options ? options : defaultoptions;
1824: /* The number of unused PETSc options */
1825: n = 0;
1826: for (i=0; i<options->N; i++) {
1827: if (!options->used[i]) n++;
1828: }
1829: if (N) { *N = n; }
1830: if (names) { PetscMalloc1(n,names); }
1831: if (values) { PetscMalloc1(n,values); }
1833: n = 0;
1834: if (names || values) {
1835: for (i=0; i<options->N; i++) {
1836: if (!options->used[i]) {
1837: if (names) (*names)[n] = options->names[i];
1838: if (values) (*values)[n] = options->values[i];
1839: n++;
1840: }
1841: }
1842: }
1843: return(0);
1844: }
1846: /*@C
1847: PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using PetscOptionsLeftGet.
1849: Not Collective
1851: Input Parameter:
1852: + options - options database, use NULL for default global database
1853: . names - names of options not used
1854: - values - values of options not used
1856: Level: advanced
1858: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet()
1859: @*/
1860: PetscErrorCode PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[])
1861: {
1868: if (N) { *N = 0; }
1869: if (names) { PetscFree(*names); }
1870: if (values) { PetscFree(*values); }
1871: return(0);
1872: }
1874: /*@C
1875: PetscOptionsMonitorDefault - Print all options set value events using the supplied PetscViewer.
1877: Logically Collective on ctx
1879: Input Parameters:
1880: + name - option name string
1881: . value - option value string
1882: - ctx - an ASCII viewer or NULL
1884: Level: intermediate
1886: Notes:
1887: If ctx=NULL, PetscPrintf() is used.
1888: The first MPI rank in the PetscViewer viewer actually prints the values, other
1889: processes may have different values set
1891: .seealso: PetscOptionsMonitorSet()
1892: @*/
1893: PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx)
1894: {
1898: if (ctx) {
1899: PetscViewer viewer = (PetscViewer)ctx;
1900: if (!value) {
1901: PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);
1902: } else if (!value[0]) {
1903: PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);
1904: } else {
1905: PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);
1906: }
1907: } else {
1908: MPI_Comm comm = PETSC_COMM_WORLD;
1909: if (!value) {
1910: PetscPrintf(comm,"Removing option: %s\n",name,value);
1911: } else if (!value[0]) {
1912: PetscPrintf(comm,"Setting option: %s (no value)\n",name);
1913: } else {
1914: PetscPrintf(comm,"Setting option: %s = %s\n",name,value);
1915: }
1916: }
1917: return(0);
1918: }
1920: /*@C
1921: PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
1922: modified the PETSc options database.
1924: Not Collective
1926: Input Parameters:
1927: + monitor - pointer to function (if this is NULL, it turns off monitoring
1928: . mctx - [optional] context for private data for the
1929: monitor routine (use NULL if no context is desired)
1930: - monitordestroy - [optional] routine that frees monitor context
1931: (may be NULL)
1933: Calling Sequence of monitor:
1934: $ monitor (const char name[], const char value[], void *mctx)
1936: + name - option name string
1937: . value - option value string
1938: - mctx - optional monitoring context, as set by PetscOptionsMonitorSet()
1940: Options Database Keys:
1941: See PetscInitialize() for options related to option database monitoring.
1943: Notes:
1944: The default is to do nothing. To print the name and value of options
1945: being inserted into the database, use PetscOptionsMonitorDefault() as the monitoring routine,
1946: with a null monitoring context.
1948: Several different monitoring routines may be set by calling
1949: PetscOptionsMonitorSet() multiple times; all will be called in the
1950: order in which they were set.
1952: Level: intermediate
1954: .seealso: PetscOptionsMonitorDefault(), PetscInitialize()
1955: @*/
1956: PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1957: {
1958: PetscOptions options = defaultoptions;
1961: if (options->monitorCancel) return(0);
1962: if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set");
1963: options->monitor[options->numbermonitors] = monitor;
1964: options->monitordestroy[options->numbermonitors] = monitordestroy;
1965: options->monitorcontext[options->numbermonitors++] = (void*)mctx;
1966: return(0);
1967: }
1969: /*
1970: PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
1971: Empty string is considered as true.
1972: */
1973: PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a)
1974: {
1975: PetscBool istrue,isfalse;
1976: size_t len;
1980: /* PetscStrlen() returns 0 for NULL or "" */
1981: PetscStrlen(value,&len);
1982: if (!len) {*a = PETSC_TRUE; return(0);}
1983: PetscStrcasecmp(value,"TRUE",&istrue);
1984: if (istrue) {*a = PETSC_TRUE; return(0);}
1985: PetscStrcasecmp(value,"YES",&istrue);
1986: if (istrue) {*a = PETSC_TRUE; return(0);}
1987: PetscStrcasecmp(value,"1",&istrue);
1988: if (istrue) {*a = PETSC_TRUE; return(0);}
1989: PetscStrcasecmp(value,"on",&istrue);
1990: if (istrue) {*a = PETSC_TRUE; return(0);}
1991: PetscStrcasecmp(value,"FALSE",&isfalse);
1992: if (isfalse) {*a = PETSC_FALSE; return(0);}
1993: PetscStrcasecmp(value,"NO",&isfalse);
1994: if (isfalse) {*a = PETSC_FALSE; return(0);}
1995: PetscStrcasecmp(value,"0",&isfalse);
1996: if (isfalse) {*a = PETSC_FALSE; return(0);}
1997: PetscStrcasecmp(value,"off",&isfalse);
1998: if (isfalse) {*a = PETSC_FALSE; return(0);}
1999: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value);
2000: }
2002: /*
2003: PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2004: */
2005: PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a)
2006: {
2008: size_t len;
2009: PetscBool decide,tdefault,mouse;
2012: PetscStrlen(name,&len);
2013: if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value");
2015: PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);
2016: if (!tdefault) {
2017: PetscStrcasecmp(name,"DEFAULT",&tdefault);
2018: }
2019: PetscStrcasecmp(name,"PETSC_DECIDE",&decide);
2020: if (!decide) {
2021: PetscStrcasecmp(name,"DECIDE",&decide);
2022: }
2023: PetscStrcasecmp(name,"mouse",&mouse);
2025: if (tdefault) *a = PETSC_DEFAULT;
2026: else if (decide) *a = PETSC_DECIDE;
2027: else if (mouse) *a = -1;
2028: else {
2029: char *endptr;
2030: long strtolval;
2032: strtolval = strtol(name,&endptr,10);
2033: if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no integer value (do not include . in it)",name);
2035: #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2036: (void) strtolval;
2037: *a = atoll(name);
2038: #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2039: (void) strtolval;
2040: *a = _atoi64(name);
2041: #else
2042: *a = (PetscInt)strtolval;
2043: #endif
2044: }
2045: return(0);
2046: }
2048: #if defined(PETSC_USE_REAL___FLOAT128)
2049: #include <quadmath.h>
2050: #endif
2052: static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr)
2053: {
2055: #if defined(PETSC_USE_REAL___FLOAT128)
2056: *a = strtoflt128(name,endptr);
2057: #else
2058: *a = (PetscReal)strtod(name,endptr);
2059: #endif
2060: return(0);
2061: }
2063: static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary)
2064: {
2065: PetscBool hasi = PETSC_FALSE;
2066: char *ptr;
2067: PetscReal strtoval;
2071: PetscStrtod(name,&strtoval,&ptr);
2072: if (ptr == name) {
2073: strtoval = 1.;
2074: hasi = PETSC_TRUE;
2075: if (name[0] == 'i') {
2076: ptr++;
2077: } else if (name[0] == '+' && name[1] == 'i') {
2078: ptr += 2;
2079: } else if (name[0] == '-' && name[1] == 'i') {
2080: strtoval = -1.;
2081: ptr += 2;
2082: }
2083: } else if (*ptr == 'i') {
2084: hasi = PETSC_TRUE;
2085: ptr++;
2086: }
2087: *endptr = ptr;
2088: *isImaginary = hasi;
2089: if (hasi) {
2090: #if !defined(PETSC_USE_COMPLEX)
2091: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name);
2092: #else
2093: *a = PetscCMPLX(0.,strtoval);
2094: #endif
2095: } else {
2096: *a = strtoval;
2097: }
2098: return(0);
2099: }
2101: /*
2102: Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2103: */
2104: PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a)
2105: {
2106: size_t len;
2107: PetscBool match;
2108: char *endptr;
2112: PetscStrlen(name,&len);
2113: if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"String of length zero has no numerical value");
2115: PetscStrcasecmp(name,"PETSC_DEFAULT",&match);
2116: if (!match) {
2117: PetscStrcasecmp(name,"DEFAULT",&match);
2118: }
2119: if (match) {*a = PETSC_DEFAULT; return(0);}
2121: PetscStrcasecmp(name,"PETSC_DECIDE",&match);
2122: if (!match) {
2123: PetscStrcasecmp(name,"DECIDE",&match);
2124: }
2125: if (match) {*a = PETSC_DECIDE; return(0);}
2127: PetscStrtod(name,a,&endptr);
2128: if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name);
2129: return(0);
2130: }
2132: PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a)
2133: {
2134: PetscBool imag1;
2135: size_t len;
2136: PetscScalar val = 0.;
2137: char *ptr = NULL;
2141: PetscStrlen(name,&len);
2142: if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value");
2143: PetscStrtoz(name,&val,&ptr,&imag1);
2144: #if defined(PETSC_USE_COMPLEX)
2145: if ((size_t) (ptr - name) < len) {
2146: PetscBool imag2;
2147: PetscScalar val2;
2149: PetscStrtoz(ptr,&val2,&ptr,&imag2);
2150: if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name);
2151: val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2));
2152: }
2153: #endif
2154: if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name);
2155: *a = val;
2156: return(0);
2157: }
2159: /*@C
2160: PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2161: option in the database.
2163: Not Collective
2165: Input Parameters:
2166: + options - options database, use NULL for default global database
2167: . pre - the string to prepend to the name or NULL
2168: - name - the option one is seeking
2170: Output Parameter:
2171: + ivalue - the logical value to return
2172: - set - PETSC_TRUE if found, else PETSC_FALSE
2174: Level: beginner
2176: Notes:
2177: TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2178: FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
2180: If the option is given, but no value is provided, then ivalue and set are both given the value PETSC_TRUE. That is -requested_bool
2181: is equivalent to -requested_bool true
2183: If the user does not supply the option at all ivalue is NOT changed. Thus
2184: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2186: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
2187: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(),
2188: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2189: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2190: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2191: PetscOptionsFList(), PetscOptionsEList()
2192: @*/
2193: PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set)
2194: {
2195: const char *value;
2196: PetscBool flag;
2202: PetscOptionsFindPair(options,pre,name,&value,&flag);
2203: if (flag) {
2204: if (set) *set = PETSC_TRUE;
2205: PetscOptionsStringToBool(value, &flag);
2206: if (ivalue) *ivalue = flag;
2207: } else {
2208: if (set) *set = PETSC_FALSE;
2209: }
2210: return(0);
2211: }
2213: /*@C
2214: PetscOptionsGetEList - Puts a list of option values that a single one may be selected from
2216: Not Collective
2218: Input Parameters:
2219: + options - options database, use NULL for default global database
2220: . pre - the string to prepend to the name or NULL
2221: . opt - option name
2222: . list - the possible choices (one of these must be selected, anything else is invalid)
2223: - ntext - number of choices
2225: Output Parameter:
2226: + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2227: - set - PETSC_TRUE if found, else PETSC_FALSE
2229: Level: intermediate
2231: Notes:
2232: If the user does not supply the option value is NOT changed. Thus
2233: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2235: See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
2237: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2238: PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2239: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2240: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2241: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2242: PetscOptionsFList(), PetscOptionsEList()
2243: @*/
2244: PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set)
2245: {
2247: size_t alen,len = 0, tlen = 0;
2248: char *svalue;
2249: PetscBool aset,flg = PETSC_FALSE;
2250: PetscInt i;
2254: for (i=0; i<ntext; i++) {
2255: PetscStrlen(list[i],&alen);
2256: if (alen > len) len = alen;
2257: tlen += len + 1;
2258: }
2259: len += 5; /* a little extra space for user mistypes */
2260: PetscMalloc1(len,&svalue);
2261: PetscOptionsGetString(options,pre,opt,svalue,len,&aset);
2262: if (aset) {
2263: PetscEListFind(ntext,list,svalue,value,&flg);
2264: if (!flg) {
2265: char *avail,*pavl;
2267: PetscMalloc1(tlen,&avail);
2268: pavl = avail;
2269: for (i=0; i<ntext; i++) {
2270: PetscStrlen(list[i],&alen);
2271: PetscStrcpy(pavl,list[i]);
2272: pavl += alen;
2273: PetscStrcpy(pavl," ");
2274: pavl += 1;
2275: }
2276: PetscStrtolower(avail);
2277: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s. Available options: %s",svalue,pre ? pre : "",opt+1,avail);
2278: }
2279: if (set) *set = PETSC_TRUE;
2280: } else if (set) *set = PETSC_FALSE;
2281: PetscFree(svalue);
2282: return(0);
2283: }
2285: /*@C
2286: PetscOptionsGetEnum - Gets the enum value for a particular option in the database.
2288: Not Collective
2290: Input Parameters:
2291: + options - options database, use NULL for default global database
2292: . pre - option prefix or NULL
2293: . opt - option name
2294: . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2295: - defaultv - the default (current) value
2297: Output Parameter:
2298: + value - the value to return
2299: - set - PETSC_TRUE if found, else PETSC_FALSE
2301: Level: beginner
2303: Notes:
2304: If the user does not supply the option value is NOT changed. Thus
2305: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2307: List is usually something like PCASMTypes or some other predefined list of enum names
2309: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2310: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2311: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
2312: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2313: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2314: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2315: PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2316: @*/
2317: PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set)
2318: {
2320: PetscInt ntext = 0,tval;
2321: PetscBool fset;
2325: while (list[ntext++]) {
2326: if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
2327: }
2328: if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
2329: ntext -= 3;
2330: PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);
2331: /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2332: if (fset) *value = (PetscEnum)tval;
2333: if (set) *set = fset;
2334: return(0);
2335: }
2337: /*@C
2338: PetscOptionsGetInt - Gets the integer value for a particular option in the database.
2340: Not Collective
2342: Input Parameters:
2343: + options - options database, use NULL for default global database
2344: . pre - the string to prepend to the name or NULL
2345: - name - the option one is seeking
2347: Output Parameter:
2348: + ivalue - the integer value to return
2349: - set - PETSC_TRUE if found, else PETSC_FALSE
2351: Level: beginner
2353: Notes:
2354: If the user does not supply the option ivalue is NOT changed. Thus
2355: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2357: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
2358: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2359: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
2360: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2361: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2362: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2363: PetscOptionsFList(), PetscOptionsEList()
2364: @*/
2365: PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set)
2366: {
2367: const char *value;
2369: PetscBool flag;
2374: PetscOptionsFindPair(options,pre,name,&value,&flag);
2375: if (flag) {
2376: if (!value) {
2377: if (set) *set = PETSC_FALSE;
2378: } else {
2379: if (set) *set = PETSC_TRUE;
2380: PetscOptionsStringToInt(value,ivalue);
2381: }
2382: } else {
2383: if (set) *set = PETSC_FALSE;
2384: }
2385: return(0);
2386: }
2388: /*@C
2389: PetscOptionsGetReal - Gets the double precision value for a particular
2390: option in the database.
2392: Not Collective
2394: Input Parameters:
2395: + options - options database, use NULL for default global database
2396: . pre - string to prepend to each name or NULL
2397: - name - the option one is seeking
2399: Output Parameter:
2400: + dvalue - the double value to return
2401: - set - PETSC_TRUE if found, PETSC_FALSE if not found
2403: Notes:
2404: If the user does not supply the option dvalue is NOT changed. Thus
2405: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2407: Level: beginner
2409: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2410: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
2411: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2412: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2413: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2414: PetscOptionsFList(), PetscOptionsEList()
2415: @*/
2416: PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set)
2417: {
2418: const char *value;
2419: PetscBool flag;
2425: PetscOptionsFindPair(options,pre,name,&value,&flag);
2426: if (flag) {
2427: if (!value) {
2428: if (set) *set = PETSC_FALSE;
2429: } else {
2430: if (set) *set = PETSC_TRUE;
2431: PetscOptionsStringToReal(value,dvalue);
2432: }
2433: } else {
2434: if (set) *set = PETSC_FALSE;
2435: }
2436: return(0);
2437: }
2439: /*@C
2440: PetscOptionsGetScalar - Gets the scalar value for a particular
2441: option in the database.
2443: Not Collective
2445: Input Parameters:
2446: + options - options database, use NULL for default global database
2447: . pre - string to prepend to each name or NULL
2448: - name - the option one is seeking
2450: Output Parameter:
2451: + dvalue - the double value to return
2452: - set - PETSC_TRUE if found, else PETSC_FALSE
2454: Level: beginner
2456: Usage:
2457: A complex number 2+3i must be specified with NO spaces
2459: Notes:
2460: If the user does not supply the option dvalue is NOT changed. Thus
2461: you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.
2463: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2464: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2465: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2466: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2467: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2468: PetscOptionsFList(), PetscOptionsEList()
2469: @*/
2470: PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set)
2471: {
2472: const char *value;
2473: PetscBool flag;
2479: PetscOptionsFindPair(options,pre,name,&value,&flag);
2480: if (flag) {
2481: if (!value) {
2482: if (set) *set = PETSC_FALSE;
2483: } else {
2484: #if !defined(PETSC_USE_COMPLEX)
2485: PetscOptionsStringToReal(value,dvalue);
2486: #else
2487: PetscOptionsStringToScalar(value,dvalue);
2488: #endif
2489: if (set) *set = PETSC_TRUE;
2490: }
2491: } else { /* flag */
2492: if (set) *set = PETSC_FALSE;
2493: }
2494: return(0);
2495: }
2497: /*@C
2498: PetscOptionsGetString - Gets the string value for a particular option in
2499: the database.
2501: Not Collective
2503: Input Parameters:
2504: + options - options database, use NULL for default global database
2505: . pre - string to prepend to name or NULL
2506: . name - the option one is seeking
2507: - len - maximum length of the string including null termination
2509: Output Parameters:
2510: + string - location to copy string
2511: - set - PETSC_TRUE if found, else PETSC_FALSE
2513: Level: beginner
2515: Fortran Note:
2516: The Fortran interface is slightly different from the C/C++
2517: interface (len is not used). Sample usage in Fortran follows
2518: .vb
2519: character *20 string
2520: PetscErrorCode ierr
2521: PetscBool set
2522: call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2523: .ve
2525: Notes:
2526: 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
2528: If the user does not use the option then the string is not changed. Thus
2529: you should ALWAYS initialize the string if you access it without first checking if the set flag is true.
2531: Note:
2532: 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).
2534: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2535: PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2536: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2537: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2538: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2539: PetscOptionsFList(), PetscOptionsEList()
2540: @*/
2541: PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set)
2542: {
2543: const char *value;
2544: PetscBool flag;
2550: PetscOptionsFindPair(options,pre,name,&value,&flag);
2551: if (!flag) {
2552: if (set) *set = PETSC_FALSE;
2553: } else {
2554: if (set) *set = PETSC_TRUE;
2555: if (value) {
2556: PetscStrncpy(string,value,len);
2557: } else {
2558: PetscArrayzero(string,len);
2559: }
2560: }
2561: return(0);
2562: }
2564: char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[])
2565: {
2566: const char *value;
2567: PetscBool flag;
2571: PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(NULL);
2572: if (flag) PetscFunctionReturn((char*)value);
2573: else PetscFunctionReturn(NULL);
2574: }
2576: /*@C
2577: PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2578: option in the database. The values must be separated with commas with
2579: no intervening spaces.
2581: Not Collective
2583: Input Parameters:
2584: + options - options database, use NULL for default global database
2585: . pre - string to prepend to each name or NULL
2586: . name - the option one is seeking
2587: - nmax - maximum number of values to retrieve
2589: Output Parameter:
2590: + dvalue - the integer values to return
2591: . nmax - actual number of values retreived
2592: - set - PETSC_TRUE if found, else PETSC_FALSE
2594: Level: beginner
2596: Notes:
2597: TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2598: FALSE, false, NO, no, and 0 all translate to PETSC_FALSE
2600: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2601: PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2602: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2603: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2604: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2605: PetscOptionsFList(), PetscOptionsEList()
2606: @*/
2607: PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set)
2608: {
2609: const char *svalue;
2610: char *value;
2612: PetscInt n = 0;
2613: PetscBool flag;
2614: PetscToken token;
2621: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2622: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2623: if (set) *set = PETSC_TRUE;
2624: PetscTokenCreate(svalue,',',&token);
2625: PetscTokenFind(token,&value);
2626: while (value && n < *nmax) {
2627: PetscOptionsStringToBool(value,dvalue);
2628: PetscTokenFind(token,&value);
2629: dvalue++;
2630: n++;
2631: }
2632: PetscTokenDestroy(&token);
2633: *nmax = n;
2634: return(0);
2635: }
2637: /*@C
2638: PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.
2640: Not Collective
2642: Input Parameters:
2643: + options - options database, use NULL for default global database
2644: . pre - option prefix or NULL
2645: . name - option name
2646: . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2647: - nmax - maximum number of values to retrieve
2649: Output Parameters:
2650: + ivalue - the enum values to return
2651: . nmax - actual number of values retreived
2652: - set - PETSC_TRUE if found, else PETSC_FALSE
2654: Level: beginner
2656: Notes:
2657: The array must be passed as a comma separated list.
2659: There must be no intervening spaces between the values.
2661: list is usually something like PCASMTypes or some other predefined list of enum names.
2663: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2664: PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2665: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(),
2666: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(),
2667: PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2668: PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2669: @*/
2670: PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set)
2671: {
2672: const char *svalue;
2673: char *value;
2674: PetscInt n = 0;
2675: PetscEnum evalue;
2676: PetscBool flag;
2677: PetscToken token;
2686: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2687: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2688: if (set) *set = PETSC_TRUE;
2689: PetscTokenCreate(svalue,',',&token);
2690: PetscTokenFind(token,&value);
2691: while (value && n < *nmax) {
2692: PetscEnumFind(list,value,&evalue,&flag);
2693: if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1);
2694: ivalue[n++] = evalue;
2695: PetscTokenFind(token,&value);
2696: }
2697: PetscTokenDestroy(&token);
2698: *nmax = n;
2699: return(0);
2700: }
2702: /*@C
2703: PetscOptionsGetIntArray - Gets an array of integer values for a particular
2704: option in the database.
2706: Not Collective
2708: Input Parameters:
2709: + options - options database, use NULL for default global database
2710: . pre - string to prepend to each name or NULL
2711: . name - the option one is seeking
2712: - nmax - maximum number of values to retrieve
2714: Output Parameter:
2715: + ivalue - the integer values to return
2716: . nmax - actual number of values retreived
2717: - set - PETSC_TRUE if found, else PETSC_FALSE
2719: Level: beginner
2721: Notes:
2722: The array can be passed as
2723: a comma separated list: 0,1,2,3,4,5,6,7
2724: a range (start-end+1): 0-8
2725: a range with given increment (start-end+1:inc): 0-7:2
2726: a combination of values and ranges separated by commas: 0,1-8,8-15:2
2728: There must be no intervening spaces between the values.
2730: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2731: PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2732: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2733: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2734: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2735: PetscOptionsFList(), PetscOptionsEList()
2736: @*/
2737: PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set)
2738: {
2739: const char *svalue;
2740: char *value;
2742: PetscInt n = 0,i,j,start,end,inc,nvalues;
2743: size_t len;
2744: PetscBool flag,foundrange;
2745: PetscToken token;
2752: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2753: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2754: if (set) *set = PETSC_TRUE;
2755: PetscTokenCreate(svalue,',',&token);
2756: PetscTokenFind(token,&value);
2757: while (value && n < *nmax) {
2758: /* look for form d-D where d and D are integers */
2759: foundrange = PETSC_FALSE;
2760: PetscStrlen(value,&len);
2761: if (value[0] == '-') i=2;
2762: else i=1;
2763: for (;i<(int)len; i++) {
2764: if (value[i] == '-') {
2765: if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
2766: value[i] = 0;
2768: PetscOptionsStringToInt(value,&start);
2769: inc = 1;
2770: j = i+1;
2771: for (;j<(int)len; j++) {
2772: if (value[j] == ':') {
2773: value[j] = 0;
2775: PetscOptionsStringToInt(value+j+1,&inc);
2776: 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);
2777: break;
2778: }
2779: }
2780: PetscOptionsStringToInt(value+i+1,&end);
2781: 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);
2782: nvalues = (end-start)/inc + (end-start)%inc;
2783: 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);
2784: for (;start<end; start+=inc) {
2785: *ivalue = start; ivalue++;n++;
2786: }
2787: foundrange = PETSC_TRUE;
2788: break;
2789: }
2790: }
2791: if (!foundrange) {
2792: PetscOptionsStringToInt(value,ivalue);
2793: ivalue++;
2794: n++;
2795: }
2796: PetscTokenFind(token,&value);
2797: }
2798: PetscTokenDestroy(&token);
2799: *nmax = n;
2800: return(0);
2801: }
2803: /*@C
2804: PetscOptionsGetRealArray - Gets an array of double precision values for a
2805: particular option in the database. The values must be separated with
2806: commas with no intervening spaces.
2808: Not Collective
2810: Input Parameters:
2811: + options - options database, use NULL for default global database
2812: . pre - string to prepend to each name or NULL
2813: . name - the option one is seeking
2814: - nmax - maximum number of values to retrieve
2816: Output Parameters:
2817: + dvalue - the double values to return
2818: . nmax - actual number of values retreived
2819: - set - PETSC_TRUE if found, else PETSC_FALSE
2821: Level: beginner
2823: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2824: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
2825: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2826: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2827: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2828: PetscOptionsFList(), PetscOptionsEList()
2829: @*/
2830: PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set)
2831: {
2832: const char *svalue;
2833: char *value;
2835: PetscInt n = 0;
2836: PetscBool flag;
2837: PetscToken token;
2844: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2845: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2846: if (set) *set = PETSC_TRUE;
2847: PetscTokenCreate(svalue,',',&token);
2848: PetscTokenFind(token,&value);
2849: while (value && n < *nmax) {
2850: PetscOptionsStringToReal(value,dvalue++);
2851: PetscTokenFind(token,&value);
2852: n++;
2853: }
2854: PetscTokenDestroy(&token);
2855: *nmax = n;
2856: return(0);
2857: }
2859: /*@C
2860: PetscOptionsGetScalarArray - Gets an array of scalars for a
2861: particular option in the database. The values must be separated with
2862: commas with no intervening spaces.
2864: Not Collective
2866: Input Parameters:
2867: + options - options database, use NULL for default global database
2868: . pre - string to prepend to each name or NULL
2869: . name - the option one is seeking
2870: - nmax - maximum number of values to retrieve
2872: Output Parameters:
2873: + dvalue - the scalar values to return
2874: . nmax - actual number of values retreived
2875: - set - PETSC_TRUE if found, else PETSC_FALSE
2877: Level: beginner
2879: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2880: PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
2881: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2882: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2883: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2884: PetscOptionsFList(), PetscOptionsEList()
2885: @*/
2886: PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],PetscInt *nmax,PetscBool *set)
2887: {
2888: const char *svalue;
2889: char *value;
2891: PetscInt n = 0;
2892: PetscBool flag;
2893: PetscToken token;
2900: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2901: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2902: if (set) *set = PETSC_TRUE;
2903: PetscTokenCreate(svalue,',',&token);
2904: PetscTokenFind(token,&value);
2905: while (value && n < *nmax) {
2906: PetscOptionsStringToScalar(value,dvalue++);
2907: PetscTokenFind(token,&value);
2908: n++;
2909: }
2910: PetscTokenDestroy(&token);
2911: *nmax = n;
2912: return(0);
2913: }
2915: /*@C
2916: PetscOptionsGetStringArray - Gets an array of string values for a particular
2917: option in the database. The values must be separated with commas with
2918: no intervening spaces.
2920: Not Collective
2922: Input Parameters:
2923: + options - options database, use NULL for default global database
2924: . pre - string to prepend to name or NULL
2925: . name - the option one is seeking
2926: - nmax - maximum number of strings
2928: Output Parameter:
2929: + strings - location to copy strings
2930: - set - PETSC_TRUE if found, else PETSC_FALSE
2932: Level: beginner
2934: Notes:
2935: The user should pass in an array of pointers to char, to hold all the
2936: strings returned by this function.
2938: The user is responsible for deallocating the strings that are
2939: returned. The Fortran interface for this routine is not supported.
2941: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2942: PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2943: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2944: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2945: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2946: PetscOptionsFList(), PetscOptionsEList()
2947: @*/
2948: PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set)
2949: {
2950: const char *svalue;
2951: char *value;
2953: PetscInt n = 0;
2954: PetscBool flag;
2955: PetscToken token;
2962: PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2963: if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2964: if (set) *set = PETSC_TRUE;
2965: PetscTokenCreate(svalue,',',&token);
2966: PetscTokenFind(token,&value);
2967: while (value && n < *nmax) {
2968: PetscStrallocpy(value,&strings[n]);
2969: PetscTokenFind(token,&value);
2970: n++;
2971: }
2972: PetscTokenDestroy(&token);
2973: *nmax = n;
2974: return(0);
2975: }
2977: /*@C
2978: PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one
2980: Prints a deprecation warning, unless an option is supplied to suppress.
2982: Logically Collective
2984: Input Parameters:
2985: + pre - string to prepend to name or NULL
2986: . oldname - the old, deprecated option
2987: . newname - the new option, or NULL if option is purely removed
2988: . version - a string describing the version of first deprecation, e.g. "3.9"
2989: - info - additional information string, or NULL.
2991: Options Database Keys:
2992: . -options_suppress_deprecated_warnings - do not print deprecation warnings
2994: Notes:
2995: Must be called between PetscOptionsBegin() (or PetscObjectOptionsBegin()) and PetscOptionsEnd().
2996: Only the proces of rank zero that owns the PetscOptionsItems are argument (managed by PetscOptionsBegin() or
2997: PetscObjectOptionsBegin() prints the information
2998: If newname is provided, the old option is replaced. Otherwise, it remains
2999: in the options database.
3000: If an option is not replaced, the info argument should be used to advise the user
3001: on how to proceed.
3002: There is a limit on the length of the warning printed, so very long strings
3003: provided as info may be truncated.
3005: Level: developer
3007: .seealso: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsScalar(), PetscOptionsBool(), PetscOptionsString(), PetscOptionsSetValue()
3009: @*/
3010: PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[])
3011: {
3012: PetscErrorCode ierr;
3013: PetscBool found,quiet;
3014: const char *value;
3015: const char * const quietopt="-options_suppress_deprecated_warnings";
3016: char msg[4096];
3017: char *prefix = NULL;
3018: PetscOptions options = NULL;
3019: MPI_Comm comm = PETSC_COMM_SELF;
3024: if (PetscOptionsObject) {
3025: prefix = PetscOptionsObject->prefix;
3026: options = PetscOptionsObject->options;
3027: comm = PetscOptionsObject->comm;
3028: }
3029: PetscOptionsFindPair(options,prefix,oldname,&value,&found);
3030: if (found) {
3031: if (newname) {
3032: if (prefix) {
3033: PetscOptionsPrefixPush(options,prefix);
3034: }
3035: PetscOptionsSetValue(options,newname,value);
3036: if (prefix) {
3037: PetscOptionsPrefixPop(options);
3038: }
3039: PetscOptionsClearValue(options,oldname);
3040: }
3041: quiet = PETSC_FALSE;
3042: PetscOptionsGetBool(options,NULL,quietopt,&quiet,NULL);
3043: if (!quiet) {
3044: PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");
3045: PetscStrcat(msg,oldname);
3046: PetscStrcat(msg," is deprecated as of version ");
3047: PetscStrcat(msg,version);
3048: PetscStrcat(msg," and will be removed in a future release.");
3049: if (newname) {
3050: PetscStrcat(msg," Please use the option ");
3051: PetscStrcat(msg,newname);
3052: PetscStrcat(msg," instead.");
3053: }
3054: if (info) {
3055: PetscStrcat(msg," ");
3056: PetscStrcat(msg,info);
3057: }
3058: PetscStrcat(msg," (Silence this warning with ");
3059: PetscStrcat(msg,quietopt);
3060: PetscStrcat(msg,")\n");
3061: PetscPrintf(comm,msg);
3062: }
3063: }
3064: return(0);
3065: }