Actual source code: options.c

  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

 22: #if defined(PETSC_HAVE_STRCASECMP)
 23: #define PetscOptNameCmp(a,b) strcasecmp(a,b)
 24: #elif defined(PETSC_HAVE_STRICMP)
 25: #define PetscOptNameCmp(a,b) stricmp(a,b)
 26: #else
 27: #define PetscOptNameCmp(a,b) Error_strcasecmp_not_found
 28: #endif

 30: #include <petsc/private/hashtable.h>

 32: /* This assumes ASCII encoding and ignores locale settings */
 33: /* Using tolower() is about 2X slower in microbenchmarks   */
 34: PETSC_STATIC_INLINE int PetscToLower(int c)
 35: {
 36:   return ((c >= 'A') & (c <= 'Z')) ? c + 'a' - 'A' : c;
 37: }

 39: /* Bob Jenkins's one at a time hash function (case-insensitive) */
 40: PETSC_STATIC_INLINE unsigned int PetscOptHash(const char key[])
 41: {
 42:   unsigned int hash = 0;
 43:   while (*key) {
 44:     hash += PetscToLower(*key++);
 45:     hash += hash << 10;
 46:     hash ^= hash >>  6;
 47:   }
 48:   hash += hash <<  3;
 49:   hash ^= hash >> 11;
 50:   hash += hash << 15;
 51:   return hash;
 52: }

 54: PETSC_STATIC_INLINE int PetscOptEqual(const char a[],const char b[])
 55: {
 56:   return !PetscOptNameCmp(a,b);
 57: }

 59: KHASH_INIT(HO, kh_cstr_t, int, 1, PetscOptHash, PetscOptEqual)

 61: /*
 62:     This table holds all the options set by the user. For simplicity, we use a static size database
 63: */
 64: #define MAXOPTNAME PETSC_MAX_OPTION_NAME
 65: #define MAXOPTIONS 512
 66: #define MAXALIASES  25
 67: #define MAXPREFIXES 25
 68: #define MAXOPTIONSMONITORS 5

 70: struct  _n_PetscOptions {
 71:   PetscOptions   previous;
 72:   int            N;                    /* number of options */
 73:   char           *names[MAXOPTIONS];   /* option names */
 74:   char           *values[MAXOPTIONS];  /* option values */
 75:   PetscBool      used[MAXOPTIONS];     /* flag option use */
 76:   PetscBool      precedentProcessed;

 78:   /* Hash table */
 79:   khash_t(HO)    *ht;

 81:   /* Prefixes */
 82:   int            prefixind;
 83:   int            prefixstack[MAXPREFIXES];
 84:   char           prefix[MAXOPTNAME];

 86:   /* Aliases */
 87:   int            Naliases;                   /* number or aliases */
 88:   char           *aliases1[MAXALIASES];      /* aliased */
 89:   char           *aliases2[MAXALIASES];      /* aliasee */

 91:   /* Help */
 92:   PetscBool      help;       /* flag whether "-help" is in the database */
 93:   PetscBool      help_intro; /* flag whether "-help intro" is in the database */

 95:   /* Monitors */
 96:   PetscBool      monitorFromOptions, monitorCancel;
 97:   PetscErrorCode (*monitor[MAXOPTIONSMONITORS])(const char[],const char[],void*); /* returns control to user after */
 98:   PetscErrorCode (*monitordestroy[MAXOPTIONSMONITORS])(void**);         /* */
 99:   void           *monitorcontext[MAXOPTIONSMONITORS];                  /* to pass arbitrary user data into monitor */
100:   PetscInt       numbermonitors;                                       /* to, for instance, detect options being set */
101: };

103: static PetscOptions defaultoptions = NULL;  /* the options database routines query this object for options */

105: /* list of options which preceed others, i.e., are processed in PetscOptionsProcessPrecedentFlags() */
106: static const char *precedentOptions[] = {"-options_monitor","-options_monitor_cancel","-help","-skip_petscrc"};
107: enum PetscPrecedentOption {PO_OPTIONS_MONITOR,PO_OPTIONS_MONITOR_CANCEL,PO_HELP,PO_SKIP_PETSCRC,PO_NUM};

109: static PetscErrorCode PetscOptionsSetValue_Private(PetscOptions,const char[],const char[],int*);

111: /*
112:     Options events monitor
113: */
114: static PetscErrorCode PetscOptionsMonitor(PetscOptions options,const char name[],const char value[])
115: {
116:   PetscInt       i;

119:   if (!PetscErrorHandlingInitialized) return 0;
121:   if (!value) value = "";
122:   if (options->monitorFromOptions) {
123:     PetscOptionsMonitorDefault(name,value,NULL);
124:   }
125:   for (i=0; i<options->numbermonitors; i++) {
126:     (*options->monitor[i])(name,value,options->monitorcontext[i]);
127:   }
128:   return(0);
129: }

131: /*@
132:    PetscOptionsCreate - Creates an empty options database.

134:    Logically collective

136:    Output Parameter:
137: .  options - Options database object

139:    Level: advanced

141:    Developer Note: We may want eventually to pass a MPI_Comm to determine the ownership of the object

143: .seealso: PetscOptionsDestroy(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsInsert(), PetscOptionsSetValue()
144: @*/
145: PetscErrorCode PetscOptionsCreate(PetscOptions *options)
146: {
147:   if (!options) return PETSC_ERR_ARG_NULL;
148:   *options = (PetscOptions)calloc(1,sizeof(**options));
149:   if (!*options) return PETSC_ERR_MEM;
150:   return 0;
151: }

153: /*@
154:     PetscOptionsDestroy - Destroys an option database.

156:     Logically collective on whatever communicator was associated with the call to PetscOptionsCreate()

158:   Input Parameter:
159: .  options - the PetscOptions object

161:    Level: advanced

163: .seealso: PetscOptionsInsert(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsInsert(), PetscOptionsSetValue()
164: @*/
165: PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
166: {

170:   if (!*options) return 0;
171:   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()");
172:   PetscOptionsClear(*options);if (ierr) return ierr;
173:   /* XXX what about monitors ? */
174:   free(*options);
175:   *options = NULL;
176:   return(0);
177: }

179: /*
180:     PetscOptionsCreateDefault - Creates the default global options database
181: */
182: PetscErrorCode PetscOptionsCreateDefault(void)
183: {

186:   if (!defaultoptions) {
187:     PetscOptionsCreate(&defaultoptions);if (ierr) return ierr;
188:   }
189:   return 0;
190: }

192: /*@
193:       PetscOptionsPush - Push a new PetscOptions object as the default provider of options
194:                          Allows using different parts of a code to use different options databases

196:   Logically Collective

198:   Input Parameter:
199: .   opt - the options obtained with PetscOptionsCreate()

201:   Notes:
202:   Use PetscOptionsPop() to return to the previous default options database

204:   The collectivity of this routine is complex; only the MPI processes that call this routine will
205:   have the affect of these options. If some processes that create objects call this routine and others do
206:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
207:   on different ranks.

209:    Level: advanced

211: .seealso: PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsInsert(), PetscOptionsSetValue(), PetscOptionsLeft()

213: @*/
214: PetscErrorCode PetscOptionsPush(PetscOptions opt)
215: {

219:   PetscOptionsCreateDefault();
220:   opt->previous  = defaultoptions;
221:   defaultoptions = opt;
222:   return(0);
223: }

225: /*@
226:       PetscOptionsPop - Pop the most recent PetscOptionsPush() to return to the previous default options

228:       Logically collective on whatever communicator was associated with the call to PetscOptionsCreate()

230:   Notes:
231:   Use PetscOptionsPop() to return to the previous default options database
232:   Allows using different parts of a code to use different options databases

234:    Level: advanced

236: .seealso: PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsInsert(), PetscOptionsSetValue(), PetscOptionsLeft()

238: @*/
239: PetscErrorCode PetscOptionsPop(void)
240: {
241:   PetscOptions current = defaultoptions;

244:   if (!defaultoptions) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing default options");
245:   if (!defaultoptions->previous) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PetscOptionsPop() called too many times");
246:   defaultoptions    = defaultoptions->previous;
247:   current->previous = NULL;
248:   return(0);
249: }

251: /*
252:     PetscOptionsDestroyDefault - Destroys the default global options database
253: */
254: PetscErrorCode PetscOptionsDestroyDefault(void)
255: {
257:   PetscOptions   tmp;

259:   if (!defaultoptions) return 0;
260:   /* Destroy any options that the user forgot to pop */
261:   while (defaultoptions->previous) {
262:     tmp = defaultoptions;
263:     PetscOptionsPop();
264:     PetscOptionsDestroy(&tmp);
265:   }
266:   PetscOptionsDestroy(&defaultoptions);if (ierr) return ierr;
267:   return 0;
268: }

270: /*@C
271:    PetscOptionsValidKey - PETSc Options database keys must begin with one or two dashes (-) followed by a letter.

273:    Not collective

275:    Input Parameter:
276: .  key - string to check if valid

278:    Output Parameter:
279: .  valid - PETSC_TRUE if a valid key

281:    Level: intermediate
282: @*/
283: PetscErrorCode PetscOptionsValidKey(const char key[],PetscBool *valid)
284: {
285:   char           *ptr;

290:   *valid = PETSC_FALSE;
291:   if (!key) return(0);
292:   if (key[0] != '-') return(0);
293:   if (key[1] == '-') key++;
294:   if (!isalpha((int)key[1])) return(0);
295:   (void) strtod(key,&ptr);
296:   if (ptr != key && !(*ptr == '_' || isalnum((int)*ptr))) return(0);
297:   *valid = PETSC_TRUE;
298:   return(0);
299: }

301: /*@C
302:    PetscOptionsInsertString - Inserts options into the database from a string

304:    Logically Collective

306:    Input Parameters:
307: +  options - options object
308: -  in_str - string that contains options separated by blanks

310:    Level: intermediate

312:   The collectivity of this routine is complex; only the MPI processes that call this routine will
313:   have the affect of these options. If some processes that create objects call this routine and others do
314:   not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
315:   on different ranks.

317:    Contributed by Boyana Norris

319: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
320:           PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
321:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
322:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
323:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
324:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsInsertFile()
325: @*/
326: PetscErrorCode PetscOptionsInsertString(PetscOptions options,const char in_str[])
327: {
328:   MPI_Comm       comm = PETSC_COMM_SELF;
330:   char           *first,*second;
331:   PetscToken     token;

334:   PetscTokenCreate(in_str,' ',&token);
335:   PetscTokenFind(token,&first);
336:   while (first) {
337:     PetscBool isfile,isfileyaml,isstringyaml,ispush,ispop,key;
338:     PetscStrcasecmp(first,"-options_file",&isfile);
339:     PetscStrcasecmp(first,"-options_file_yaml",&isfileyaml);
340:     PetscStrcasecmp(first,"-options_string_yaml",&isstringyaml);
341:     PetscStrcasecmp(first,"-prefix_push",&ispush);
342:     PetscStrcasecmp(first,"-prefix_pop",&ispop);
343:     PetscOptionsValidKey(first,&key);
344:     if (!key) {
345:       PetscTokenFind(token,&first);
346:     } else if (isfile) {
347:       PetscTokenFind(token,&second);
348:       PetscOptionsInsertFile(comm,options,second,PETSC_TRUE);
349:       PetscTokenFind(token,&first);
350:     } else if (isfileyaml) {
351:       PetscTokenFind(token,&second);
352:       PetscOptionsInsertFileYAML(comm,options,second,PETSC_TRUE);
353:       PetscTokenFind(token,&first);
354:     } else if (isstringyaml) {
355:       PetscTokenFind(token,&second);
356:       PetscOptionsInsertStringYAML(options,second);
357:       PetscTokenFind(token,&first);
358:     } else if (ispush) {
359:       PetscTokenFind(token,&second);
360:       PetscOptionsPrefixPush(options,second);
361:       PetscTokenFind(token,&first);
362:     } else if (ispop) {
363:       PetscOptionsPrefixPop(options);
364:       PetscTokenFind(token,&first);
365:     } else {
366:       PetscTokenFind(token,&second);
367:       PetscOptionsValidKey(second,&key);
368:       if (!key) {
369:         PetscOptionsSetValue(options,first,second);
370:         PetscTokenFind(token,&first);
371:       } else {
372:         PetscOptionsSetValue(options,first,NULL);
373:         first = second;
374:       }
375:     }
376:   }
377:   PetscTokenDestroy(&token);
378:   return(0);
379: }

381: /*
382:     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
383: */
384: static char *Petscgetline(FILE * f)
385: {
386:   size_t size  = 0;
387:   size_t len   = 0;
388:   size_t last  = 0;
389:   char   *buf  = NULL;

391:   if (feof(f)) return NULL;
392:   do {
393:     size += 1024; /* BUFSIZ is defined as "the optimal read size for this platform" */
394:     buf   = (char*)realloc((void*)buf,size); /* realloc(NULL,n) is the same as malloc(n) */
395:     /* Actually do the read. Note that fgets puts a terminal '\0' on the
396:     end of the string, so we make sure we overwrite this */
397:     if (!fgets(buf+len,1024,f)) buf[len]=0;
398:     PetscStrlen(buf,&len);
399:     last = len - 1;
400:   } while (!feof(f) && buf[last] != '\n' && buf[last] != '\r');
401:   if (len) return buf;
402:   free(buf);
403:   return NULL;
404: }

406: static PetscErrorCode PetscOptionsFilename(MPI_Comm comm,const char file[],char filename[PETSC_MAX_PATH_LEN],PetscBool *yaml)
407: {
408:   char           fname[PETSC_MAX_PATH_LEN+8],path[PETSC_MAX_PATH_LEN+8],*tail;

412:   *yaml = PETSC_FALSE;
413:   PetscStrreplace(comm,file,fname,sizeof(fname));
414:   PetscFixFilename(fname,path);
415:   PetscStrendswith(path,":yaml",yaml);
416:   if (*yaml) {
417:     PetscStrrchr(path,':',&tail);
418:     tail[-1] = 0; /* remove ":yaml" suffix from path */
419:   }
420:   PetscStrncpy(filename,path,PETSC_MAX_PATH_LEN);
421:   /* check for standard YAML and JSON filename extensions */
422:   if (!*yaml) {PetscStrendswith(filename,".yaml",yaml);}
423:   if (!*yaml) {PetscStrendswith(filename,".yml", yaml);}
424:   if (!*yaml) {PetscStrendswith(filename,".json",yaml);}
425:   if (!*yaml) { /* check file contents */
426:     PetscMPIInt rank;
427:     MPI_Comm_rank(comm,&rank);
428:     if (rank == 0) {
429:       FILE *fh = fopen(filename,"r");
430:       if (fh) {
431:         char buf[6] = "";
432:         if (fread(buf,1,6,fh) > 0) {
433:           PetscStrncmp(buf,"%YAML ",6,yaml);  /* check for '%YAML' tag */
434:           if (!*yaml) {PetscStrncmp(buf,"---",3,yaml);}  /* check for document start */
435:         }
436:         (void)fclose(fh);
437:       }
438:     }
439:     MPI_Bcast(yaml,1,MPIU_BOOL,0,comm);
440:   }
441:   return(0);
442: }

444: static PetscErrorCode PetscOptionsInsertFilePetsc(MPI_Comm comm,PetscOptions options,const char file[],PetscBool require)
445: {
446:   char           *string,*vstring = NULL,*astring = NULL,*packed = NULL;
447:   char           *tokens[4];
449:   size_t         i,len,bytes;
450:   FILE           *fd;
451:   PetscToken     token=NULL;
452:   int            err;
453:   char           *cmatch;
454:   const char     cmt='#';
455:   PetscInt       line=1;
456:   PetscMPIInt    rank,cnt=0,acnt=0,counts[2];
457:   PetscBool      isdir,alias=PETSC_FALSE,valid;

460:   PetscMemzero(tokens,sizeof(tokens));
461:   MPI_Comm_rank(comm,&rank);
462:   if (rank == 0) {
463:     char fpath[PETSC_MAX_PATH_LEN];
464:     char fname[PETSC_MAX_PATH_LEN];

466:     PetscStrreplace(PETSC_COMM_SELF,file,fpath,sizeof(fpath));
467:     PetscFixFilename(fpath,fname);

469:     fd   = fopen(fname,"r");
470:     PetscTestDirectory(fname,'r',&isdir);
471:     if (isdir && require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified options file %s is a directory",fname);
472:     if (fd && !isdir) {
473:       PetscSegBuffer vseg,aseg;
474:       PetscSegBufferCreate(1,4000,&vseg);
475:       PetscSegBufferCreate(1,2000,&aseg);

477:       /* the following line will not work when opening initial files (like .petscrc) since info is not yet set */
478:       PetscInfo1(NULL,"Opened options file %s\n",file);

480:       while ((string = Petscgetline(fd))) {
481:         /* eliminate comments from each line */
482:         PetscStrchr(string,cmt,&cmatch);
483:         if (cmatch) *cmatch = 0;
484:         PetscStrlen(string,&len);
485:         /* replace tabs, ^M, \n with " " */
486:         for (i=0; i<len; i++) {
487:           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') {
488:             string[i] = ' ';
489:           }
490:         }
491:         PetscTokenCreate(string,' ',&token);
492:         PetscTokenFind(token,&tokens[0]);
493:         if (!tokens[0]) {
494:           goto destroy;
495:         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
496:           PetscTokenFind(token,&tokens[0]);
497:         }
498:         for (i=1; i<4; i++) {
499:           PetscTokenFind(token,&tokens[i]);
500:         }
501:         if (!tokens[0]) {
502:           goto destroy;
503:         } else if (tokens[0][0] == '-') {
504:           PetscOptionsValidKey(tokens[0],&valid);
505:           if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: invalid option %s",fname,line,tokens[0]);
506:           PetscStrlen(tokens[0],&len);
507:           PetscSegBufferGet(vseg,len+1,&vstring);
508:           PetscArraycpy(vstring,tokens[0],len);
509:           vstring[len] = ' ';
510:           if (tokens[1]) {
511:             PetscOptionsValidKey(tokens[1],&valid);
512:             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]);
513:             PetscStrlen(tokens[1],&len);
514:             PetscSegBufferGet(vseg,len+3,&vstring);
515:             vstring[0] = '"';
516:             PetscArraycpy(vstring+1,tokens[1],len);
517:             vstring[len+1] = '"';
518:             vstring[len+2] = ' ';
519:           }
520:         } else {
521:           PetscStrcasecmp(tokens[0],"alias",&alias);
522:           if (alias) {
523:             PetscOptionsValidKey(tokens[1],&valid);
524:             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]);
525:             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]);
526:             PetscOptionsValidKey(tokens[2],&valid);
527:             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]);
528:             PetscStrlen(tokens[1],&len);
529:             PetscSegBufferGet(aseg,len+1,&astring);
530:             PetscArraycpy(astring,tokens[1],len);
531:             astring[len] = ' ';

533:             PetscStrlen(tokens[2],&len);
534:             PetscSegBufferGet(aseg,len+1,&astring);
535:             PetscArraycpy(astring,tokens[2],len);
536:             astring[len] = ' ';
537:           } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown first token in options file %s line %D: %s",fname,line,tokens[0]);
538:         }
539:         {
540:           const char *extraToken = alias ? tokens[3] : tokens[2];
541:           if (extraToken) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: extra token %s",fname,line,extraToken);
542:         }
543: destroy:
544:         free(string);
545:         PetscTokenDestroy(&token);
546:         alias = PETSC_FALSE;
547:         line++;
548:       }
549:       err = fclose(fd);
550:       if (err) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file %s",fname);
551:       PetscSegBufferGetSize(aseg,&bytes); /* size without null termination */
552:       PetscMPIIntCast(bytes,&acnt);
553:       PetscSegBufferGet(aseg,1,&astring);
554:       astring[0] = 0;
555:       PetscSegBufferGetSize(vseg,&bytes); /* size without null termination */
556:       PetscMPIIntCast(bytes,&cnt);
557:       PetscSegBufferGet(vseg,1,&vstring);
558:       vstring[0] = 0;
559:       PetscMalloc1(2+acnt+cnt,&packed);
560:       PetscSegBufferExtractTo(aseg,packed);
561:       PetscSegBufferExtractTo(vseg,packed+acnt+1);
562:       PetscSegBufferDestroy(&aseg);
563:       PetscSegBufferDestroy(&vseg);
564:     } else if (require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Unable to open options file %s",fname);
565:   }

567:   counts[0] = acnt;
568:   counts[1] = cnt;
569:   err = MPI_Bcast(counts,2,MPI_INT,0,comm);
570:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in first MPI collective call, could be caused by using an incorrect mpiexec or a network problem, it can be caused by having VPN running: see https://petsc.org/release/faq/");
571:   acnt = counts[0];
572:   cnt = counts[1];
573:   if (rank) {
574:     PetscMalloc1(2+acnt+cnt,&packed);
575:   }
576:   if (acnt || cnt) {
577:     MPI_Bcast(packed,2+acnt+cnt,MPI_CHAR,0,comm);
578:     astring = packed;
579:     vstring = packed + acnt + 1;
580:   }

582:   if (acnt) {
583:     PetscTokenCreate(astring,' ',&token);
584:     PetscTokenFind(token,&tokens[0]);
585:     while (tokens[0]) {
586:       PetscTokenFind(token,&tokens[1]);
587:       PetscOptionsSetAlias(options,tokens[0],tokens[1]);
588:       PetscTokenFind(token,&tokens[0]);
589:     }
590:     PetscTokenDestroy(&token);
591:   }

593:   if (cnt) {
594:     PetscOptionsInsertString(options,vstring);
595:   }
596:   PetscFree(packed);
597:   return(0);
598: }

600: /*@C
601:      PetscOptionsInsertFile - Inserts options into the database from a file.

603:      Collective

605:   Input Parameters:
606: +   comm - the processes that will share the options (usually PETSC_COMM_WORLD)
607: .   options - options database, use NULL for default global database
608: .   file - name of file,
609:            ".yml" and ".yaml" filename extensions are inserted as YAML options,
610:            append ":yaml" to filename to force YAML options.
611: -   require - if PETSC_TRUE will generate an error if the file does not exist

613:   Notes:
614:    Use  # for lines that are comments and which should be ignored.
615:    Usually, instead of using this command, one should list the file name in the call to PetscInitialize(), this insures that certain options
616:    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
617:    calls to XXXSetFromOptions() it should not be used for options listed under PetscInitialize().
618:    The collectivity of this routine is complex; only the MPI processes in comm will
619:    have the affect of these options. If some processes that create objects call this routine and others do
620:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
621:    on different ranks.

623:   Level: developer

625: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
626:           PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
627:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
628:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
629:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
630:           PetscOptionsFList(), PetscOptionsEList()

632: @*/
633: PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm,PetscOptions options,const char file[],PetscBool require)
634: {
635:   char           filename[PETSC_MAX_PATH_LEN];
636:   PetscBool      yaml;

640:   PetscOptionsFilename(comm,file,filename,&yaml);
641:   if (yaml) {
642:     PetscOptionsInsertFileYAML(comm,options,filename,require);
643:   } else {
644:     PetscOptionsInsertFilePetsc(comm,options,filename,require);
645:   }
646:   return(0);
647: }

649: /*@C
650:    PetscOptionsInsertArgs - Inserts options into the database from a array of strings

652:    Logically Collective

654:    Input Parameters:
655: +  options - options object
656: .  argc - the array lenght
657: -  args - the string array

659:    Level: intermediate

661: .seealso: PetscOptions, PetscOptionsInsertString(), PetscOptionsInsertFile()
662: @*/
663: PetscErrorCode PetscOptionsInsertArgs(PetscOptions options,int argc,char *args[])
664: {
665:   MPI_Comm       comm = PETSC_COMM_WORLD;
667:   int            left          = PetscMax(argc,0);
668:   char           *const *eargs = args;

671:   while (left) {
672:     PetscBool isfile,isfileyaml,isstringyaml,ispush,ispop,key;
673:     PetscStrcasecmp(eargs[0],"-options_file",&isfile);
674:     PetscStrcasecmp(eargs[0],"-options_file_yaml",&isfileyaml);
675:     PetscStrcasecmp(eargs[0],"-options_string_yaml",&isstringyaml);
676:     PetscStrcasecmp(eargs[0],"-prefix_push",&ispush);
677:     PetscStrcasecmp(eargs[0],"-prefix_pop",&ispop);
678:     PetscOptionsValidKey(eargs[0],&key);
679:     if (!key) {
680:       eargs++; left--;
681:     } else if (isfile) {
682:       if (left <= 1 || eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file filename option");
683:       PetscOptionsInsertFile(comm,options,eargs[1],PETSC_TRUE);
684:       eargs += 2; left -= 2;
685:     } else if (isfileyaml) {
686:       if (left <= 1 || eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file_yaml filename option");
687:       PetscOptionsInsertFileYAML(comm,options,eargs[1],PETSC_TRUE);
688:       eargs += 2; left -= 2;
689:     } else if (isstringyaml) {
690:       if (left <= 1 || eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing string for -options_string_yaml string option");
691:       PetscOptionsInsertStringYAML(options,eargs[1]);
692:       eargs += 2; left -= 2;
693:     } else if (ispush) {
694:       if (left <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option");
695:       if (eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option (prefixes cannot start with '-')");
696:       PetscOptionsPrefixPush(options,eargs[1]);
697:       eargs += 2; left -= 2;
698:     } else if (ispop) {
699:       PetscOptionsPrefixPop(options);
700:       eargs++; left--;
701:     } else {
702:       PetscBool nextiskey = PETSC_FALSE;
703:       if (left >= 2) {PetscOptionsValidKey(eargs[1],&nextiskey);}
704:       if (left < 2 || nextiskey) {
705:         PetscOptionsSetValue(options,eargs[0],NULL);
706:         eargs++; left--;
707:       } else {
708:         PetscOptionsSetValue(options,eargs[0],eargs[1]);
709:         eargs += 2; left -= 2;
710:       }
711:     }
712:   }
713:   return(0);
714: }

716: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsStringToBoolIfSet_Private(enum PetscPrecedentOption opt,const char *val[],PetscBool set[],PetscBool *flg)
717: {

721:   if (set[opt]) {
722:     PetscOptionsStringToBool(val[opt],flg);
723:   } else *flg = PETSC_FALSE;
724:   return(0);
725: }

727: /* Process options with absolute precedence */
728: static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options,int argc,char *args[],PetscBool *skip_petscrc,PetscBool *skip_petscrc_set)
729: {
730:   const char* const *opt = precedentOptions;
731:   const size_t      n = PO_NUM;
732:   size_t            o;
733:   int               a;
734:   const char        **val;
735:   PetscBool         *set;
736:   PetscErrorCode    ierr;

739:   PetscCalloc2(n,&val,n,&set);

741:   /* Look for options possibly set using PetscOptionsSetValue beforehand */
742:   for (o=0; o<n; o++) {
743:     PetscOptionsFindPair(options,NULL,opt[o],&val[o],&set[o]);
744:   }

746:   /* Loop through all args to collect last occurring value of each option */
747:   for (a=1; a<argc; a++) {
748:     PetscBool valid, eq;

750:     PetscOptionsValidKey(args[a],&valid);
751:     if (!valid) continue;
752:     for (o=0; o<n; o++) {
753:       PetscStrcasecmp(args[a],opt[o],&eq);
754:       if (eq) {
755:         set[o] = PETSC_TRUE;
756:         if (a == argc-1 || !args[a+1] || !args[a+1][0] || args[a+1][0] == '-') val[o] = NULL;
757:         else val[o] = args[a+1];
758:         break;
759:       }
760:     }
761:   }

763:   /* Process flags */
764:   PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro);
765:   if (options->help_intro) options->help = PETSC_TRUE;
766:   else {PetscOptionsStringToBoolIfSet_Private(PO_HELP,            val,set,&options->help);}
767:   PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL,val,set,&options->monitorCancel);
768:   PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR,       val,set,&options->monitorFromOptions);
769:   PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC,          val,set,skip_petscrc);
770:   *skip_petscrc_set = set[PO_SKIP_PETSCRC];

772:   /* Store precedent options in database and mark them as used */
773:   for (o=0; o<n; o++) {
774:     if (set[o]) {
775:       PetscOptionsSetValue_Private(options,opt[o],val[o],&a);
776:       options->used[a] = PETSC_TRUE;
777:     }
778:   }

780:   PetscFree2(val,set);
781:   options->precedentProcessed = PETSC_TRUE;
782:   return(0);
783: }

785: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options,const char name[],PetscBool *flg)
786: {
787:   int i;

790:   *flg = PETSC_FALSE;
791:   if (options->precedentProcessed) {
792:     for (i=0; i<PO_NUM; i++) {
793:       if (!PetscOptNameCmp(precedentOptions[i],name)) {
794:         /* check if precedent option has been set already */
795:         PetscOptionsFindPair(options,NULL,name,NULL,flg);if (ierr) return ierr;
796:         if (*flg) break;
797:       }
798:     }
799:   }
800:   return 0;
801: }

803: /*@C
804:    PetscOptionsInsert - Inserts into the options database from the command line,
805:                         the environmental variable and a file.

807:    Collective on PETSC_COMM_WORLD

809:    Input Parameters:
810: +  options - options database or NULL for the default global database
811: .  argc - count of number of command line arguments
812: .  args - the command line arguments
813: -  file - [optional] PETSc database file, append ":yaml" to filename to specify YAML options format.
814:           Use NULL or empty string to not check for code specific file.
815:           Also checks ~/.petscrc, .petscrc and petscrc.
816:           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.

818:    Note:
819:    Since PetscOptionsInsert() is automatically called by PetscInitialize(),
820:    the user does not typically need to call this routine. PetscOptionsInsert()
821:    can be called several times, adding additional entries into the database.

823:    Options Database Keys:
824: +   -options_file <filename> - read options from a file
825: -   -options_file_yaml <filename> - read options from a YAML file

827:    See PetscInitialize() for options related to option database monitoring.

829:    Level: advanced

831: .seealso: PetscOptionsDestroy(), PetscOptionsView(), PetscOptionsInsertString(), PetscOptionsInsertFile(),
832:           PetscInitialize()
833: @*/
834: PetscErrorCode PetscOptionsInsert(PetscOptions options,int *argc,char ***args,const char file[])
835: {
836:   MPI_Comm       comm = PETSC_COMM_WORLD;
838:   PetscMPIInt    rank;
839:   PetscBool      hasArgs = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
840:   PetscBool      skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;

843:   if (hasArgs && !(args && *args)) SETERRQ(comm,PETSC_ERR_ARG_NULL,"*argc > 1 but *args not given");
844:   MPI_Comm_rank(comm,&rank);

846:   if (!options) {
847:     PetscOptionsCreateDefault();
848:     options = defaultoptions;
849:   }
850:   if (hasArgs) {
851:     /* process options with absolute precedence */
852:     PetscOptionsProcessPrecedentFlags(options,*argc,*args,&skipPetscrc,&skipPetscrcSet);
853:   }
854:   if (file && file[0]) {
855:     PetscOptionsInsertFile(comm,options,file,PETSC_TRUE);
856:     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
857:     if (!skipPetscrcSet) {PetscOptionsGetBool(options,NULL,"-skip_petscrc",&skipPetscrc,NULL);}
858:   }
859:   if (!skipPetscrc) {
860:     char filename[PETSC_MAX_PATH_LEN];
861:     PetscGetHomeDirectory(filename,sizeof(filename));
862:     MPI_Bcast(filename,(int)sizeof(filename),MPI_CHAR,0,comm);
863:     if (filename[0]) {PetscStrcat(filename,"/.petscrc");}
864:     PetscOptionsInsertFile(comm,options,filename,PETSC_FALSE);
865:     PetscOptionsInsertFile(comm,options,".petscrc",PETSC_FALSE);
866:     PetscOptionsInsertFile(comm,options,"petscrc",PETSC_FALSE);
867:   }

869:   /* insert environment options */
870:   {
871:     char   *eoptions = NULL;
872:     size_t len       = 0;
873:     if (rank == 0) {
874:       eoptions = (char*)getenv("PETSC_OPTIONS");
875:       PetscStrlen(eoptions,&len);
876:     }
877:     MPI_Bcast(&len,1,MPIU_SIZE_T,0,comm);
878:     if (len) {
879:       if (rank) {PetscMalloc1(len+1,&eoptions);}
880:       MPI_Bcast(eoptions,len,MPI_CHAR,0,comm);
881:       if (rank) eoptions[len] = 0;
882:       PetscOptionsInsertString(options,eoptions);
883:       if (rank) {PetscFree(eoptions);}
884:     }
885:   }

887:   /* insert YAML environment options */
888:   {
889:     char   *eoptions = NULL;
890:     size_t len       = 0;
891:     if (rank == 0) {
892:       eoptions = (char*)getenv("PETSC_OPTIONS_YAML");
893:       PetscStrlen(eoptions,&len);
894:     }
895:     MPI_Bcast(&len,1,MPIU_SIZE_T,0,comm);
896:     if (len) {
897:       if (rank) {PetscMalloc1(len+1,&eoptions);}
898:       MPI_Bcast(eoptions,len,MPI_CHAR,0,comm);
899:       if (rank) eoptions[len] = 0;
900:       PetscOptionsInsertStringYAML(options,eoptions);
901:       if (rank) {PetscFree(eoptions);}
902:     }
903:   }

905:   /* insert command line options here because they take precedence over arguments in petscrc/environment */
906:   if (hasArgs) {PetscOptionsInsertArgs(options,*argc-1,*args+1);}
907:   return(0);
908: }

910: /*@C
911:    PetscOptionsView - Prints the options that have been loaded. This is
912:    useful for debugging purposes.

914:    Logically Collective on PetscViewer

916:    Input Parameters:
917: +  options - options database, use NULL for default global database
918: -  viewer - must be an PETSCVIEWERASCII viewer

920:    Options Database Key:
921: .  -options_view - Activates PetscOptionsView() within PetscFinalize()

923:    Notes:
924:    Only the rank zero process of MPI_Comm used to create view prints the option values. Other processes
925:    may have different values but they are not printed.

927:    Level: advanced

929: .seealso: PetscOptionsAllUsed()
930: @*/
931: PetscErrorCode PetscOptionsView(PetscOptions options,PetscViewer viewer)
932: {
934:   PetscInt       i;
935:   PetscBool      isascii;

939:   options = options ? options : defaultoptions;
940:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
941:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
942:   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Only supports ASCII viewer");

944:   if (!options->N) {
945:     PetscViewerASCIIPrintf(viewer,"#No PETSc Option Table entries\n");
946:     return(0);
947:   }

949:   PetscViewerASCIIPrintf(viewer,"#PETSc Option Table entries:\n");
950:   for (i=0; i<options->N; i++) {
951:     if (options->values[i]) {
952:       PetscViewerASCIIPrintf(viewer,"-%s %s\n",options->names[i],options->values[i]);
953:     } else {
954:       PetscViewerASCIIPrintf(viewer,"-%s\n",options->names[i]);
955:     }
956:   }
957:   PetscViewerASCIIPrintf(viewer,"#End of PETSc Option Table entries\n");
958:   return(0);
959: }

961: /*
962:    Called by error handlers to print options used in run
963: */
964: PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
965: {
966:   PetscInt     i;
967:   PetscOptions options = defaultoptions;

970:   if (options->N) {
971:     (*PetscErrorPrintf)("PETSc Option Table entries:\n");
972:   } else {
973:     (*PetscErrorPrintf)("No PETSc Option Table entries\n");
974:   }
975:   for (i=0; i<options->N; i++) {
976:     if (options->values[i]) {
977:       (*PetscErrorPrintf)("-%s %s\n",options->names[i],options->values[i]);
978:     } else {
979:       (*PetscErrorPrintf)("-%s\n",options->names[i]);
980:     }
981:   }
982:   return(0);
983: }

985: /*@C
986:    PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow.

988:    Logically Collective

990:    Input Parameters:
991: +  options - options database, or NULL for the default global database
992: -  prefix - The string to append to the existing prefix

994:    Options Database Keys:
995: +   -prefix_push <some_prefix_> - push the given prefix
996: -   -prefix_pop - pop the last prefix

998:    Notes:
999:    It is common to use this in conjunction with -options_file as in

1001: $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop

1003:    where the files no longer require all options to be prefixed with -system2_.

1005:    The collectivity of this routine is complex; only the MPI processes that call this routine will
1006:    have the affect of these options. If some processes that create objects call this routine and others do
1007:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1008:    on different ranks.

1010: Level: advanced

1012: .seealso: PetscOptionsPrefixPop(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
1013: @*/
1014: PetscErrorCode PetscOptionsPrefixPush(PetscOptions options,const char prefix[])
1015: {
1017:   size_t         n;
1018:   PetscInt       start;
1019:   char           key[MAXOPTNAME+1];
1020:   PetscBool      valid;

1024:   options = options ? options : defaultoptions;
1025:   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);
1026:   key[0] = '-'; /* keys must start with '-' */
1027:   PetscStrncpy(key+1,prefix,sizeof(key)-1);
1028:   PetscOptionsValidKey(key,&valid);
1029:   if (!valid && options->prefixind > 0 && isdigit((int)prefix[0])) valid = PETSC_TRUE; /* If the prefix stack is not empty, make numbers a valid prefix */
1030:   if (!valid) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Given prefix \"%s\" not valid (the first character must be a letter%s, do not include leading '-')",prefix,options->prefixind?" or digit":"");
1031:   start = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
1032:   PetscStrlen(prefix,&n);
1033:   if (n+1 > sizeof(options->prefix)-start) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Maximum prefix length %d exceeded",sizeof(options->prefix));
1034:   PetscArraycpy(options->prefix+start,prefix,n+1);
1035:   options->prefixstack[options->prefixind++] = start+n;
1036:   return(0);
1037: }

1039: /*@C
1040:    PetscOptionsPrefixPop - Remove the latest options prefix, see PetscOptionsPrefixPush() for details

1042:    Logically Collective on the MPI_Comm that called PetscOptionsPrefixPush()

1044:   Input Parameters:
1045: .  options - options database, or NULL for the default global database

1047:    Level: advanced

1049: .seealso: PetscOptionsPrefixPush(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
1050: @*/
1051: PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1052: {
1053:   PetscInt offset;

1056:   options = options ? options : defaultoptions;
1057:   if (options->prefixind < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"More prefixes popped than pushed");
1058:   options->prefixind--;
1059:   offset = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
1060:   options->prefix[offset] = 0;
1061:   return(0);
1062: }

1064: /*@C
1065:     PetscOptionsClear - Removes all options form the database leaving it empty.

1067:     Logically Collective

1069:   Input Parameters:
1070: .  options - options database, use NULL for the default global database

1072:    The collectivity of this routine is complex; only the MPI processes that call this routine will
1073:    have the affect of these options. If some processes that create objects call this routine and others do
1074:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1075:    on different ranks.

1077:    Level: developer

1079: .seealso: PetscOptionsInsert()
1080: @*/
1081: PetscErrorCode PetscOptionsClear(PetscOptions options)
1082: {
1083:   PetscInt i;

1085:   options = options ? options : defaultoptions;
1086:   if (!options) return 0;

1088:   for (i=0; i<options->N; i++) {
1089:     if (options->names[i])  free(options->names[i]);
1090:     if (options->values[i]) free(options->values[i]);
1091:   }
1092:   options->N = 0;

1094:   for (i=0; i<options->Naliases; i++) {
1095:     free(options->aliases1[i]);
1096:     free(options->aliases2[i]);
1097:   }
1098:   options->Naliases = 0;

1100:   /* destroy hash table */
1101:   kh_destroy(HO,options->ht);
1102:   options->ht = NULL;

1104:   options->prefixind = 0;
1105:   options->prefix[0] = 0;
1106:   options->help      = PETSC_FALSE;
1107:   return 0;
1108: }

1110: /*@C
1111:    PetscOptionsSetAlias - Makes a key and alias for another key

1113:    Logically Collective

1115:    Input Parameters:
1116: +  options - options database, or NULL for default global database
1117: .  newname - the alias
1118: -  oldname - the name that alias will refer to

1120:    Level: advanced

1122:    The collectivity of this routine is complex; only the MPI processes that call this routine will
1123:    have the affect of these options. If some processes that create objects call this routine and others do
1124:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1125:    on different ranks.

1127: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1128:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
1129:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1130:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1131:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1132:           PetscOptionsFList(), PetscOptionsEList()
1133: @*/
1134: PetscErrorCode PetscOptionsSetAlias(PetscOptions options,const char newname[],const char oldname[])
1135: {
1136:   PetscInt       n;
1137:   size_t         len;
1138:   PetscBool      valid;

1144:   options = options ? options : defaultoptions;
1145:   PetscOptionsValidKey(newname,&valid);
1146:   if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliased option %s",newname);
1147:   PetscOptionsValidKey(oldname,&valid);
1148:   if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliasee option %s",oldname);

1150:   n = options->Naliases;
1151:   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);

1153:   newname++; oldname++;
1154:   PetscStrlen(newname,&len);
1155:   options->aliases1[n] = (char*)malloc((len+1)*sizeof(char));
1156:   PetscStrcpy(options->aliases1[n],newname);
1157:   PetscStrlen(oldname,&len);
1158:   options->aliases2[n] = (char*)malloc((len+1)*sizeof(char));
1159:   PetscStrcpy(options->aliases2[n],oldname);
1160:   options->Naliases++;
1161:   return(0);
1162: }

1164: /*@C
1165:    PetscOptionsSetValue - Sets an option name-value pair in the options
1166:    database, overriding whatever is already present.

1168:    Logically Collective

1170:    Input Parameters:
1171: +  options - options database, use NULL for the default global database
1172: .  name - name of option, this SHOULD have the - prepended
1173: -  value - the option value (not used for all options, so can be NULL)

1175:    Level: intermediate

1177:    Note:
1178:    This function can be called BEFORE PetscInitialize()

1180:    The collectivity of this routine is complex; only the MPI processes that call this routine will
1181:    have the affect of these options. If some processes that create objects call this routine and others do
1182:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1183:    on different ranks.

1185:    Developers Note: Uses malloc() directly because PETSc may not be initialized yet.

1187: .seealso: PetscOptionsInsert(), PetscOptionsClearValue()
1188: @*/
1189: PetscErrorCode PetscOptionsSetValue(PetscOptions options,const char name[],const char value[])
1190: {
1191:   return PetscOptionsSetValue_Private(options,name,value,NULL);
1192: }

1194: static PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options,const char name[],const char value[],int *pos)
1195: {
1196:   size_t         len;
1197:   int            N,n,i;
1198:   char           **names;
1199:   char           fullname[MAXOPTNAME] = "";
1200:   PetscBool      flg;

1203:   if (!options) {
1204:     PetscOptionsCreateDefault();if (ierr) return ierr;
1205:     options = defaultoptions;
1206:   }

1208:   if (name[0] != '-') return PETSC_ERR_ARG_OUTOFRANGE;

1210:   PetscOptionsSkipPrecedent(options,name,&flg);if (ierr) return ierr;
1211:   if (flg) return 0;

1213:   name++; /* skip starting dash */

1215:   if (options->prefixind > 0) {
1216:     strncpy(fullname,options->prefix,sizeof(fullname));
1217:     fullname[sizeof(fullname)-1] = 0;
1218:     strncat(fullname,name,sizeof(fullname)-strlen(fullname)-1);
1219:     fullname[sizeof(fullname)-1] = 0;
1220:     name = fullname;
1221:   }

1223:   /* check against aliases */
1224:   N = options->Naliases;
1225:   for (i=0; i<N; i++) {
1226:     int result = PetscOptNameCmp(options->aliases1[i],name);
1227:     if (!result) { name = options->aliases2[i]; break; }
1228:   }

1230:   /* slow search */
1231:   N = n = options->N;
1232:   names = options->names;
1233:   for (i=0; i<N; i++) {
1234:     int result = PetscOptNameCmp(names[i],name);
1235:     if (!result) {
1236:       n = i; goto setvalue;
1237:     } else if (result > 0) {
1238:       n = i; break;
1239:     }
1240:   }
1241:   if (N >= MAXOPTIONS) return PETSC_ERR_MEM;
1242:   /* shift remaining values up 1 */
1243:   for (i=N; i>n; i--) {
1244:     options->names[i]  = options->names[i-1];
1245:     options->values[i] = options->values[i-1];
1246:     options->used[i]   = options->used[i-1];
1247:   }
1248:   options->names[n]  = NULL;
1249:   options->values[n] = NULL;
1250:   options->used[n]   = PETSC_FALSE;
1251:   options->N++;

1253:   /* destroy hash table */
1254:   kh_destroy(HO,options->ht);
1255:   options->ht = NULL;

1257:   /* set new name */
1258:   len = strlen(name);
1259:   options->names[n] = (char*)malloc((len+1)*sizeof(char));
1260:   if (!options->names[n]) return PETSC_ERR_MEM;
1261:   strcpy(options->names[n],name);

1263: setvalue:
1264:   /* set new value */
1265:   if (options->values[n]) free(options->values[n]);
1266:   len = value ? strlen(value) : 0;
1267:   if (len) {
1268:     options->values[n] = (char*)malloc((len+1)*sizeof(char));
1269:     if (!options->values[n]) return PETSC_ERR_MEM;
1270:     strcpy(options->values[n],value);
1271:   } else {
1272:     options->values[n] = NULL;
1273:   }

1275:   /* handle -help so that it can be set from anywhere */
1276:   if (!PetscOptNameCmp(name,"help")) {
1277:     options->help = PETSC_TRUE;
1278:     options->help_intro = (value && !PetscOptNameCmp(value,"intro")) ? PETSC_TRUE : PETSC_FALSE;
1279:     options->used[n] = PETSC_TRUE;
1280:   }

1282:   if (PetscErrorHandlingInitialized) {
1283:     PetscOptionsMonitor(options,name,value);
1284:   }
1285:   if (pos) *pos = n;
1286:   return 0;
1287: }

1289: /*@C
1290:    PetscOptionsClearValue - Clears an option name-value pair in the options
1291:    database, overriding whatever is already present.

1293:    Logically Collective

1295:    Input Parameters:
1296: +  options - options database, use NULL for the default global database
1297: -  name - name of option, this SHOULD have the - prepended

1299:    Level: intermediate

1301:    The collectivity of this routine is complex; only the MPI processes that call this routine will
1302:    have the affect of these options. If some processes that create objects call this routine and others do
1303:    not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options
1304:    on different ranks.

1306: .seealso: PetscOptionsInsert()
1307: @*/
1308: PetscErrorCode PetscOptionsClearValue(PetscOptions options,const char name[])
1309: {
1310:   int            N,n,i;
1311:   char           **names;

1315:   options = options ? options : defaultoptions;
1316:   if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);
1317:   if (!PetscOptNameCmp(name,"-help")) options->help = options->help_intro = PETSC_FALSE;

1319:   name++; /* skip starting dash */

1321:   /* slow search */
1322:   N = n = options->N;
1323:   names = options->names;
1324:   for (i=0; i<N; i++) {
1325:     int result = PetscOptNameCmp(names[i],name);
1326:     if (!result) {
1327:       n = i; break;
1328:     } else if (result > 0) {
1329:       n = N; break;
1330:     }
1331:   }
1332:   if (n == N) return(0); /* it was not present */

1334:   /* remove name and value */
1335:   if (options->names[n])  free(options->names[n]);
1336:   if (options->values[n]) free(options->values[n]);
1337:   /* shift remaining values down 1 */
1338:   for (i=n; i<N-1; i++) {
1339:     options->names[i]  = options->names[i+1];
1340:     options->values[i] = options->values[i+1];
1341:     options->used[i]   = options->used[i+1];
1342:   }
1343:   options->N--;

1345:   /* destroy hash table */
1346:   kh_destroy(HO,options->ht);
1347:   options->ht = NULL;

1349:   PetscOptionsMonitor(options,name,NULL);
1350:   return(0);
1351: }

1353: /*@C
1354:    PetscOptionsFindPair - Gets an option name-value pair from the options database.

1356:    Not Collective

1358:    Input Parameters:
1359: +  options - options database, use NULL for the default global database
1360: .  pre - the string to prepend to the name or NULL, this SHOULD NOT have the "-" prepended
1361: -  name - name of option, this SHOULD have the "-" prepended

1363:    Output Parameters:
1364: +  value - the option value (optional, not used for all options)
1365: -  set - whether the option is set (optional)

1367:    Notes:
1368:    Each process may find different values or no value depending on how options were inserted into the database

1370:    Level: developer

1372: .seealso: PetscOptionsSetValue(), PetscOptionsClearValue()
1373: @*/
1374: PetscErrorCode PetscOptionsFindPair(PetscOptions options,const char pre[],const char name[],const char *value[],PetscBool *set)
1375: {
1376:   char           buf[MAXOPTNAME];
1377:   PetscBool      usehashtable = PETSC_TRUE;
1378:   PetscBool      matchnumbers = PETSC_TRUE;

1382:   options = options ? options : defaultoptions;
1383:   if (pre && PetscUnlikely(pre[0] == '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre);
1384:   if (PetscUnlikely(name[0] != '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);

1386:   name++; /* skip starting dash */

1388:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1389:   if (pre && pre[0]) {
1390:     char *ptr = buf;
1391:     if (name[0] == '-') { *ptr++ = '-';  name++; }
1392:     PetscStrncpy(ptr,pre,buf+sizeof(buf)-ptr);
1393:     PetscStrlcat(buf,name,sizeof(buf));
1394:     name = buf;
1395:   }

1397:   if (PetscDefined(USE_DEBUG)) {
1398:     PetscBool valid;
1399:     char      key[MAXOPTNAME+1] = "-";
1400:     PetscStrncpy(key+1,name,sizeof(key)-1);
1401:     PetscOptionsValidKey(key,&valid);
1402:     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1403:   }

1405:   if (!options->ht && usehashtable) {
1406:     int i,ret;
1407:     khiter_t it;
1408:     khash_t(HO) *ht;
1409:     ht = kh_init(HO);
1410:     if (PetscUnlikely(!ht)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1411:     ret = kh_resize(HO,ht,options->N*2); /* twice the required size to reduce risk of collisions */
1412:     if (PetscUnlikely(ret)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1413:     for (i=0; i<options->N; i++) {
1414:       it = kh_put(HO,ht,options->names[i],&ret);
1415:       if (PetscUnlikely(ret != 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1416:       kh_val(ht,it) = i;
1417:     }
1418:     options->ht = ht;
1419:   }

1421:   if (usehashtable)
1422:   { /* fast search */
1423:     khash_t(HO) *ht = options->ht;
1424:     khiter_t it = kh_get(HO,ht,name);
1425:     if (it != kh_end(ht)) {
1426:       int i = kh_val(ht,it);
1427:       options->used[i]  = PETSC_TRUE;
1428:       if (value) *value = options->values[i];
1429:       if (set)   *set   = PETSC_TRUE;
1430:       return(0);
1431:     }
1432:   } else
1433:   { /* slow search */
1434:     int i, N = options->N;
1435:     for (i=0; i<N; i++) {
1436:       int result = PetscOptNameCmp(options->names[i],name);
1437:       if (!result) {
1438:         options->used[i]  = PETSC_TRUE;
1439:         if (value) *value = options->values[i];
1440:         if (set)   *set   = PETSC_TRUE;
1441:         return(0);
1442:       } else if (result > 0) {
1443:         break;
1444:       }
1445:     }
1446:   }

1448:   /*
1449:    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
1450:    Maybe this special lookup mode should be enabled on request with a push/pop API.
1451:    The feature of matching _%d_ used sparingly in the codebase.
1452:    */
1453:   if (matchnumbers) {
1454:     int i,j,cnt = 0,locs[16],loce[16];
1455:     /* determine the location and number of all _%d_ in the key */
1456:     for (i=0; name[i]; i++) {
1457:       if (name[i] == '_') {
1458:         for (j=i+1; name[j]; j++) {
1459:           if (name[j] >= '0' && name[j] <= '9') continue;
1460:           if (name[j] == '_' && j > i+1) { /* found a number */
1461:             locs[cnt]   = i+1;
1462:             loce[cnt++] = j+1;
1463:           }
1464:           i = j-1;
1465:           break;
1466:         }
1467:       }
1468:     }
1469:     for (i=0; i<cnt; i++) {
1470:       PetscBool found;
1471:       char      opt[MAXOPTNAME+1] = "-", tmp[MAXOPTNAME];
1472:       PetscStrncpy(tmp,name,PetscMin((size_t)(locs[i]+1),sizeof(tmp)));
1473:       PetscStrlcat(opt,tmp,sizeof(opt));
1474:       PetscStrlcat(opt,name+loce[i],sizeof(opt));
1475:       PetscOptionsFindPair(options,NULL,opt,value,&found);
1476:       if (found) {if (set) *set = PETSC_TRUE; return(0);}
1477:     }
1478:   }

1480:   if (set) *set = PETSC_FALSE;
1481:   return(0);
1482: }

1484: /* Check whether any option begins with pre+name */
1485: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[],const char *value[],PetscBool *set)
1486: {
1487:   char           buf[MAXOPTNAME];
1488:   int            numCnt = 0, locs[16],loce[16];

1492:   options = options ? options : defaultoptions;
1493:   if (pre && pre[0] == '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre);
1494:   if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name);

1496:   name++; /* skip starting dash */

1498:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1499:   if (pre && pre[0]) {
1500:     char *ptr = buf;
1501:     if (name[0] == '-') { *ptr++ = '-';  name++; }
1502:     PetscStrncpy(ptr,pre,sizeof(buf)+(size_t)(ptr-buf));
1503:     PetscStrlcat(buf,name,sizeof(buf));
1504:     name = buf;
1505:   }

1507:   if (PetscDefined(USE_DEBUG)) {
1508:     PetscBool valid;
1509:     char      key[MAXOPTNAME+1] = "-";
1510:     PetscStrncpy(key+1,name,sizeof(key)-1);
1511:     PetscOptionsValidKey(key,&valid);
1512:     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1513:   }

1515:   /* determine the location and number of all _%d_ in the key */
1516:   {
1517:     int i,j;
1518:     for (i=0; name[i]; i++) {
1519:       if (name[i] == '_') {
1520:         for (j=i+1; name[j]; j++) {
1521:           if (name[j] >= '0' && name[j] <= '9') continue;
1522:           if (name[j] == '_' && j > i+1) { /* found a number */
1523:             locs[numCnt]   = i+1;
1524:             loce[numCnt++] = j+1;
1525:           }
1526:           i = j-1;
1527:           break;
1528:         }
1529:       }
1530:     }
1531:   }

1533:   { /* slow search */
1534:     int       c, i;
1535:     size_t    len;
1536:     PetscBool match;

1538:     for (c = -1; c < numCnt; ++c) {
1539:       char opt[MAXOPTNAME+1] = "", tmp[MAXOPTNAME];

1541:       if (c < 0) {
1542:         PetscStrcpy(opt,name);
1543:       } else {
1544:         PetscStrncpy(tmp,name,PetscMin((size_t)(locs[c]+1),sizeof(tmp)));
1545:         PetscStrlcat(opt,tmp,sizeof(opt));
1546:         PetscStrlcat(opt,name+loce[c],sizeof(opt));
1547:       }
1548:       PetscStrlen(opt,&len);
1549:       for (i=0; i<options->N; i++) {
1550:         PetscStrncmp(options->names[i],opt,len,&match);
1551:         if (match) {
1552:           options->used[i]  = PETSC_TRUE;
1553:           if (value) *value = options->values[i];
1554:           if (set)   *set   = PETSC_TRUE;
1555:           return(0);
1556:         }
1557:       }
1558:     }
1559:   }

1561:   if (set) *set = PETSC_FALSE;
1562:   return(0);
1563: }

1565: /*@C
1566:    PetscOptionsReject - Generates an error if a certain option is given.

1568:    Not Collective

1570:    Input Parameters:
1571: +  options - options database, use NULL for default global database
1572: .  pre - the option prefix (may be NULL)
1573: .  name - the option name one is seeking
1574: -  mess - error message (may be NULL)

1576:    Level: advanced

1578: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1579:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1580:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1581:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1582:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1583:           PetscOptionsFList(), PetscOptionsEList()
1584: @*/
1585: PetscErrorCode PetscOptionsReject(PetscOptions options,const char pre[],const char name[],const char mess[])
1586: {
1588:   PetscBool      flag = PETSC_FALSE;

1591:   PetscOptionsHasName(options,pre,name,&flag);
1592:   if (flag) {
1593:     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);
1594:     else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s",pre?pre:"",name+1);
1595:   }
1596:   return(0);
1597: }

1599: /*@C
1600:    PetscOptionsHasHelp - Determines whether the "-help" option is in the database.

1602:    Not Collective

1604:    Input Parameters:
1605: .  options - options database, use NULL for default global database

1607:    Output Parameters:
1608: .  set - PETSC_TRUE if found else PETSC_FALSE.

1610:    Level: advanced

1612: .seealso: PetscOptionsHasName()
1613: @*/
1614: PetscErrorCode PetscOptionsHasHelp(PetscOptions options,PetscBool *set)
1615: {
1618:   options = options ? options : defaultoptions;
1619:   *set = options->help;
1620:   return(0);
1621: }

1623: PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options,PetscBool *set)
1624: {
1627:   options = options ? options : defaultoptions;
1628:   *set = options->help_intro;
1629:   return(0);
1630: }

1632: /*@C
1633:    PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or Boolean, even
1634:                       if its value is set to false.

1636:    Not Collective

1638:    Input Parameters:
1639: +  options - options database, use NULL for default global database
1640: .  pre - string to prepend to the name or NULL
1641: -  name - the option one is seeking

1643:    Output Parameters:
1644: .  set - PETSC_TRUE if found else PETSC_FALSE.

1646:    Level: beginner

1648:    Notes:
1649:    In many cases you probably want to use PetscOptionsGetBool() instead of calling this, to allowing toggling values.

1651: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1652:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1653:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1654:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1655:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1656:           PetscOptionsFList(), PetscOptionsEList()
1657: @*/
1658: PetscErrorCode PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool *set)
1659: {
1660:   const char     *value;
1662:   PetscBool      flag;

1665:   PetscOptionsFindPair(options,pre,name,&value,&flag);
1666:   if (set) *set = flag;
1667:   return(0);
1668: }

1670: /*@C
1671:    PetscOptionsGetAll - Lists all the options the program was run with in a single string.

1673:    Not Collective

1675:    Input Parameter:
1676: .  options - the options database, use NULL for the default global database

1678:    Output Parameter:
1679: .  copts - pointer where string pointer is stored

1681:    Notes:
1682:     The array and each entry in the array should be freed with PetscFree()
1683:     Each process may have different values depending on how the options were inserted into the database

1685:    Level: advanced

1687: .seealso: PetscOptionsAllUsed(), PetscOptionsView(), PetscOptionsPush(), PetscOptionsPop()
1688: @*/
1689: PetscErrorCode PetscOptionsGetAll(PetscOptions options,char *copts[])
1690: {
1692:   PetscInt       i;
1693:   size_t         len = 1,lent = 0;
1694:   char           *coptions = NULL;

1698:   options = options ? options : defaultoptions;
1699:   /* count the length of the required string */
1700:   for (i=0; i<options->N; i++) {
1701:     PetscStrlen(options->names[i],&lent);
1702:     len += 2 + lent;
1703:     if (options->values[i]) {
1704:       PetscStrlen(options->values[i],&lent);
1705:       len += 1 + lent;
1706:     }
1707:   }
1708:   PetscMalloc1(len,&coptions);
1709:   coptions[0] = 0;
1710:   for (i=0; i<options->N; i++) {
1711:     PetscStrcat(coptions,"-");
1712:     PetscStrcat(coptions,options->names[i]);
1713:     PetscStrcat(coptions," ");
1714:     if (options->values[i]) {
1715:       PetscStrcat(coptions,options->values[i]);
1716:       PetscStrcat(coptions," ");
1717:     }
1718:   }
1719:   *copts = coptions;
1720:   return(0);
1721: }

1723: /*@C
1724:    PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database

1726:    Not Collective

1728:    Input Parameters:
1729: +  options - options database, use NULL for default global database
1730: -  name - string name of option

1732:    Output Parameter:
1733: .  used - PETSC_TRUE if the option was used, otherwise false, including if option was not found in options database

1735:    Level: advanced

1737:    Notes:
1738:    The value returned may be different on each process and depends on which options have been processed
1739:    on the given process

1741: .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed()
1742: @*/
1743: PetscErrorCode PetscOptionsUsed(PetscOptions options,const char *name,PetscBool *used)
1744: {
1745:   PetscInt       i;

1751:   options = options ? options : defaultoptions;
1752:   *used = PETSC_FALSE;
1753:   for (i=0; i<options->N; i++) {
1754:     PetscStrcasecmp(options->names[i],name,used);
1755:     if (*used) {
1756:       *used = options->used[i];
1757:       break;
1758:     }
1759:   }
1760:   return(0);
1761: }

1763: /*@
1764:    PetscOptionsAllUsed - Returns a count of the number of options in the
1765:    database that have never been selected.

1767:    Not Collective

1769:    Input Parameter:
1770: .  options - options database, use NULL for default global database

1772:    Output Parameter:
1773: .  N - count of options not used

1775:    Level: advanced

1777:    Notes:
1778:    The value returned may be different on each process and depends on which options have been processed
1779:    on the given process

1781: .seealso: PetscOptionsView()
1782: @*/
1783: PetscErrorCode PetscOptionsAllUsed(PetscOptions options,PetscInt *N)
1784: {
1785:   PetscInt     i,n = 0;

1789:   options = options ? options : defaultoptions;
1790:   for (i=0; i<options->N; i++) {
1791:     if (!options->used[i]) n++;
1792:   }
1793:   *N = n;
1794:   return(0);
1795: }

1797: /*@
1798:    PetscOptionsLeft - Prints to screen any options that were set and never used.

1800:    Not Collective

1802:    Input Parameter:
1803: .  options - options database; use NULL for default global database

1805:    Options Database Key:
1806: .  -options_left - activates PetscOptionsAllUsed() within PetscFinalize()

1808:    Notes:
1809:       This is rarely used directly, it is called by PetscFinalize() in debug more or if -options_left
1810:       is passed otherwise to help users determine possible mistakes in their usage of options. This
1811:       only prints values on process zero of PETSC_COMM_WORLD. Other processes depending the objects
1812:       used may have different options that are left unused.

1814:    Level: advanced

1816: .seealso: PetscOptionsAllUsed()
1817: @*/
1818: PetscErrorCode PetscOptionsLeft(PetscOptions options)
1819: {
1821:   PetscInt       i;
1822:   PetscInt       cnt = 0;
1823:   PetscOptions   toptions;

1826:   toptions = options ? options : defaultoptions;
1827:   for (i=0; i<toptions->N; i++) {
1828:     if (!toptions->used[i]) {
1829:       if (toptions->values[i]) {
1830:         PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",toptions->names[i],toptions->values[i]);
1831:       } else {
1832:         PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",toptions->names[i]);
1833:       }
1834:     }
1835:   }
1836:   if (!options) {
1837:     toptions = defaultoptions;
1838:     while (toptions->previous) {
1839:       cnt++;
1840:       toptions = toptions->previous;
1841:     }
1842:     if (cnt) {
1843:       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);
1844:     }
1845:   }
1846:   return(0);
1847: }

1849: /*@C
1850:    PetscOptionsLeftGet - Returns all options that were set and never used.

1852:    Not Collective

1854:    Input Parameter:
1855: .  options - options database, use NULL for default global database

1857:    Output Parameters:
1858: +  N - count of options not used
1859: .  names - names of options not used
1860: -  values - values of options not used

1862:    Level: advanced

1864:    Notes:
1865:    Users should call PetscOptionsLeftRestore() to free the memory allocated in this routine
1866:    Notes: The value returned may be different on each process and depends on which options have been processed
1867:    on the given process

1869: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft()
1870: @*/
1871: PetscErrorCode PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[])
1872: {
1874:   PetscInt       i,n;

1880:   options = options ? options : defaultoptions;

1882:   /* The number of unused PETSc options */
1883:   n = 0;
1884:   for (i=0; i<options->N; i++) {
1885:     if (!options->used[i]) n++;
1886:   }
1887:   if (N) { *N = n; }
1888:   if (names)  { PetscMalloc1(n,names); }
1889:   if (values) { PetscMalloc1(n,values); }

1891:   n = 0;
1892:   if (names || values) {
1893:     for (i=0; i<options->N; i++) {
1894:       if (!options->used[i]) {
1895:         if (names)  (*names)[n]  = options->names[i];
1896:         if (values) (*values)[n] = options->values[i];
1897:         n++;
1898:       }
1899:     }
1900:   }
1901:   return(0);
1902: }

1904: /*@C
1905:    PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using PetscOptionsLeftGet.

1907:    Not Collective

1909:    Input Parameters:
1910: +  options - options database, use NULL for default global database
1911: .  names - names of options not used
1912: -  values - values of options not used

1914:    Level: advanced

1916: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet()
1917: @*/
1918: PetscErrorCode PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[])
1919: {

1926:   if (N) { *N = 0; }
1927:   if (names)  { PetscFree(*names); }
1928:   if (values) { PetscFree(*values); }
1929:   return(0);
1930: }

1932: /*@C
1933:    PetscOptionsMonitorDefault - Print all options set value events using the supplied PetscViewer.

1935:    Logically Collective on ctx

1937:    Input Parameters:
1938: +  name  - option name string
1939: .  value - option value string
1940: -  ctx - an ASCII viewer or NULL

1942:    Level: intermediate

1944:    Notes:
1945:      If ctx=NULL, PetscPrintf() is used.
1946:      The first MPI rank in the PetscViewer viewer actually prints the values, other
1947:      processes may have different values set

1949: .seealso: PetscOptionsMonitorSet()
1950: @*/
1951: PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx)
1952: {

1956:   if (ctx) {
1957:     PetscViewer viewer = (PetscViewer)ctx;
1958:     if (!value) {
1959:       PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);
1960:     } else if (!value[0]) {
1961:       PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);
1962:     } else {
1963:       PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);
1964:     }
1965:   } else {
1966:     MPI_Comm comm = PETSC_COMM_WORLD;
1967:     if (!value) {
1968:       PetscPrintf(comm,"Removing option: %s\n",name,value);
1969:     } else if (!value[0]) {
1970:       PetscPrintf(comm,"Setting option: %s (no value)\n",name);
1971:     } else {
1972:       PetscPrintf(comm,"Setting option: %s = %s\n",name,value);
1973:     }
1974:   }
1975:   return(0);
1976: }

1978: /*@C
1979:    PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that
1980:    modified the PETSc options database.

1982:    Not Collective

1984:    Input Parameters:
1985: +  monitor - pointer to function (if this is NULL, it turns off monitoring
1986: .  mctx    - [optional] context for private data for the
1987:              monitor routine (use NULL if no context is desired)
1988: -  monitordestroy - [optional] routine that frees monitor context
1989:           (may be NULL)

1991:    Calling Sequence of monitor:
1992: $     monitor (const char name[], const char value[], void *mctx)

1994: +  name - option name string
1995: .  value - option value string
1996: -  mctx  - optional monitoring context, as set by PetscOptionsMonitorSet()

1998:    Options Database Keys:
1999:    See PetscInitialize() for options related to option database monitoring.

2001:    Notes:
2002:    The default is to do nothing.  To print the name and value of options
2003:    being inserted into the database, use PetscOptionsMonitorDefault() as the monitoring routine,
2004:    with a null monitoring context.

2006:    Several different monitoring routines may be set by calling
2007:    PetscOptionsMonitorSet() multiple times; all will be called in the
2008:    order in which they were set.

2010:    Level: intermediate

2012: .seealso: PetscOptionsMonitorDefault(), PetscInitialize()
2013: @*/
2014: PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
2015: {
2016:   PetscOptions options = defaultoptions;

2019:   if (options->monitorCancel) return(0);
2020:   if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set");
2021:   options->monitor[options->numbermonitors]          = monitor;
2022:   options->monitordestroy[options->numbermonitors]   = monitordestroy;
2023:   options->monitorcontext[options->numbermonitors++] = (void*)mctx;
2024:   return(0);
2025: }

2027: /*
2028:    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
2029:      Empty string is considered as true.
2030: */
2031: PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a)
2032: {
2033:   PetscBool      istrue,isfalse;
2034:   size_t         len;

2038:   /* PetscStrlen() returns 0 for NULL or "" */
2039:   PetscStrlen(value,&len);
2040:   if (!len)  {*a = PETSC_TRUE; return(0);}
2041:   PetscStrcasecmp(value,"TRUE",&istrue);
2042:   if (istrue) {*a = PETSC_TRUE; return(0);}
2043:   PetscStrcasecmp(value,"YES",&istrue);
2044:   if (istrue) {*a = PETSC_TRUE; return(0);}
2045:   PetscStrcasecmp(value,"1",&istrue);
2046:   if (istrue) {*a = PETSC_TRUE; return(0);}
2047:   PetscStrcasecmp(value,"on",&istrue);
2048:   if (istrue) {*a = PETSC_TRUE; return(0);}
2049:   PetscStrcasecmp(value,"FALSE",&isfalse);
2050:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2051:   PetscStrcasecmp(value,"NO",&isfalse);
2052:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2053:   PetscStrcasecmp(value,"0",&isfalse);
2054:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2055:   PetscStrcasecmp(value,"off",&isfalse);
2056:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2057:   SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value);
2058: }

2060: /*
2061:    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2062: */
2063: PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a)
2064: {
2066:   size_t         len;
2067:   PetscBool      decide,tdefault,mouse;

2070:   PetscStrlen(name,&len);
2071:   if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value");

2073:   PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);
2074:   if (!tdefault) {
2075:     PetscStrcasecmp(name,"DEFAULT",&tdefault);
2076:   }
2077:   PetscStrcasecmp(name,"PETSC_DECIDE",&decide);
2078:   if (!decide) {
2079:     PetscStrcasecmp(name,"DECIDE",&decide);
2080:   }
2081:   PetscStrcasecmp(name,"mouse",&mouse);

2083:   if (tdefault)    *a = PETSC_DEFAULT;
2084:   else if (decide) *a = PETSC_DECIDE;
2085:   else if (mouse)  *a = -1;
2086:   else {
2087:     char *endptr;
2088:     long strtolval;

2090:     strtolval = strtol(name,&endptr,10);
2091:     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);

2093: #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2094:     (void) strtolval;
2095:     *a = atoll(name);
2096: #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2097:     (void) strtolval;
2098:     *a = _atoi64(name);
2099: #else
2100:     *a = (PetscInt)strtolval;
2101: #endif
2102:   }
2103:   return(0);
2104: }

2106: #if defined(PETSC_USE_REAL___FLOAT128)
2107: #include <quadmath.h>
2108: #endif

2110: static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr)
2111: {
2113: #if defined(PETSC_USE_REAL___FLOAT128)
2114:   *a = strtoflt128(name,endptr);
2115: #else
2116:   *a = (PetscReal)strtod(name,endptr);
2117: #endif
2118:   return(0);
2119: }

2121: static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary)
2122: {
2123:   PetscBool      hasi = PETSC_FALSE;
2124:   char           *ptr;
2125:   PetscReal      strtoval;

2129:   PetscStrtod(name,&strtoval,&ptr);
2130:   if (ptr == name) {
2131:     strtoval = 1.;
2132:     hasi = PETSC_TRUE;
2133:     if (name[0] == 'i') {
2134:       ptr++;
2135:     } else if (name[0] == '+' && name[1] == 'i') {
2136:       ptr += 2;
2137:     } else if (name[0] == '-' && name[1] == 'i') {
2138:       strtoval = -1.;
2139:       ptr += 2;
2140:     }
2141:   } else if (*ptr == 'i') {
2142:     hasi = PETSC_TRUE;
2143:     ptr++;
2144:   }
2145:   *endptr = ptr;
2146:   *isImaginary = hasi;
2147:   if (hasi) {
2148: #if !defined(PETSC_USE_COMPLEX)
2149:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name);
2150: #else
2151:     *a = PetscCMPLX(0.,strtoval);
2152: #endif
2153:   } else {
2154:     *a = strtoval;
2155:   }
2156:   return(0);
2157: }

2159: /*
2160:    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2161: */
2162: PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a)
2163: {
2164:   size_t         len;
2165:   PetscBool      match;
2166:   char           *endptr;

2170:   PetscStrlen(name,&len);
2171:   if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"String of length zero has no numerical value");

2173:   PetscStrcasecmp(name,"PETSC_DEFAULT",&match);
2174:   if (!match) {
2175:     PetscStrcasecmp(name,"DEFAULT",&match);
2176:   }
2177:   if (match) {*a = PETSC_DEFAULT; return(0);}

2179:   PetscStrcasecmp(name,"PETSC_DECIDE",&match);
2180:   if (!match) {
2181:     PetscStrcasecmp(name,"DECIDE",&match);
2182:   }
2183:   if (match) {*a = PETSC_DECIDE; return(0);}

2185:   PetscStrtod(name,a,&endptr);
2186:   if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name);
2187:   return(0);
2188: }

2190: PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a)
2191: {
2192:   PetscBool      imag1;
2193:   size_t         len;
2194:   PetscScalar    val = 0.;
2195:   char           *ptr = NULL;

2199:   PetscStrlen(name,&len);
2200:   if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value");
2201:   PetscStrtoz(name,&val,&ptr,&imag1);
2202: #if defined(PETSC_USE_COMPLEX)
2203:   if ((size_t) (ptr - name) < len) {
2204:     PetscBool   imag2;
2205:     PetscScalar val2;

2207:     PetscStrtoz(ptr,&val2,&ptr,&imag2);
2208:     if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name);
2209:     val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2));
2210:   }
2211: #endif
2212:   if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name);
2213:   *a = val;
2214:   return(0);
2215: }

2217: /*@C
2218:    PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2219:             option in the database.

2221:    Not Collective

2223:    Input Parameters:
2224: +  options - options database, use NULL for default global database
2225: .  pre - the string to prepend to the name or NULL
2226: -  name - the option one is seeking

2228:    Output Parameters:
2229: +  ivalue - the logical value to return
2230: -  set - PETSC_TRUE  if found, else PETSC_FALSE

2232:    Level: beginner

2234:    Notes:
2235:        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2236:        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE

2238:       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
2239:      is equivalent to -requested_bool true

2241:        If the user does not supply the option at all ivalue is NOT changed. Thus
2242:      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2244: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
2245:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(),
2246:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2247:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2248:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2249:           PetscOptionsFList(), PetscOptionsEList()
2250: @*/
2251: PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set)
2252: {
2253:   const char     *value;
2254:   PetscBool      flag;

2260:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2261:   if (flag) {
2262:     if (set) *set = PETSC_TRUE;
2263:     PetscOptionsStringToBool(value, &flag);
2264:     if (ivalue) *ivalue = flag;
2265:   } else {
2266:     if (set) *set = PETSC_FALSE;
2267:   }
2268:   return(0);
2269: }

2271: /*@C
2272:    PetscOptionsGetEList - Puts a list of option values that a single one may be selected from

2274:    Not Collective

2276:    Input Parameters:
2277: +  options - options database, use NULL for default global database
2278: .  pre - the string to prepend to the name or NULL
2279: .  opt - option name
2280: .  list - the possible choices (one of these must be selected, anything else is invalid)
2281: -  ntext - number of choices

2283:    Output Parameters:
2284: +  value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2285: -  set - PETSC_TRUE if found, else PETSC_FALSE

2287:    Level: intermediate

2289:    Notes:
2290:     If the user does not supply the option value is NOT changed. Thus
2291:      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2293:    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()

2295: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2296:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2297:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2298:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2299:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2300:           PetscOptionsFList(), PetscOptionsEList()
2301: @*/
2302: PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set)
2303: {
2305:   size_t         alen,len = 0, tlen = 0;
2306:   char           *svalue;
2307:   PetscBool      aset,flg = PETSC_FALSE;
2308:   PetscInt       i;

2312:   for (i=0; i<ntext; i++) {
2313:     PetscStrlen(list[i],&alen);
2314:     if (alen > len) len = alen;
2315:     tlen += len + 1;
2316:   }
2317:   len += 5; /* a little extra space for user mistypes */
2318:   PetscMalloc1(len,&svalue);
2319:   PetscOptionsGetString(options,pre,opt,svalue,len,&aset);
2320:   if (aset) {
2321:     PetscEListFind(ntext,list,svalue,value,&flg);
2322:     if (!flg) {
2323:       char *avail,*pavl;

2325:       PetscMalloc1(tlen,&avail);
2326:       pavl = avail;
2327:       for (i=0; i<ntext; i++) {
2328:         PetscStrlen(list[i],&alen);
2329:         PetscStrcpy(pavl,list[i]);
2330:         pavl += alen;
2331:         PetscStrcpy(pavl," ");
2332:         pavl += 1;
2333:       }
2334:       PetscStrtolower(avail);
2335:       SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s. Available options: %s",svalue,pre ? pre : "",opt+1,avail);
2336:     }
2337:     if (set) *set = PETSC_TRUE;
2338:   } else if (set) *set = PETSC_FALSE;
2339:   PetscFree(svalue);
2340:   return(0);
2341: }

2343: /*@C
2344:    PetscOptionsGetEnum - Gets the enum value for a particular option in the database.

2346:    Not Collective

2348:    Input Parameters:
2349: +  options - options database, use NULL for default global database
2350: .  pre - option prefix or NULL
2351: .  opt - option name
2352: -  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2354:    Output Parameters:
2355: +  value - the  value to return
2356: -  set - PETSC_TRUE if found, else PETSC_FALSE

2358:    Level: beginner

2360:    Notes:
2361:     If the user does not supply the option value is NOT changed. Thus
2362:      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2364:           List is usually something like PCASMTypes or some other predefined list of enum names

2366: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2367:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2368:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
2369:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2370:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2371:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2372:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2373: @*/
2374: PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set)
2375: {
2377:   PetscInt       ntext = 0,tval;
2378:   PetscBool      fset;

2382:   while (list[ntext++]) {
2383:     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
2384:   }
2385:   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
2386:   ntext -= 3;
2387:   PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);
2388:   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2389:   if (fset) *value = (PetscEnum)tval;
2390:   if (set) *set = fset;
2391:   return(0);
2392: }

2394: /*@C
2395:    PetscOptionsGetInt - Gets the integer value for a particular option in the database.

2397:    Not Collective

2399:    Input Parameters:
2400: +  options - options database, use NULL for default global database
2401: .  pre - the string to prepend to the name or NULL
2402: -  name - the option one is seeking

2404:    Output Parameters:
2405: +  ivalue - the integer value to return
2406: -  set - PETSC_TRUE if found, else PETSC_FALSE

2408:    Level: beginner

2410:    Notes:
2411:    If the user does not supply the option ivalue is NOT changed. Thus
2412:    you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2414: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
2415:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2416:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
2417:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2418:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2419:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2420:           PetscOptionsFList(), PetscOptionsEList()
2421: @*/
2422: PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set)
2423: {
2424:   const char     *value;
2426:   PetscBool      flag;

2431:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2432:   if (flag) {
2433:     if (!value) {
2434:       if (set) *set = PETSC_FALSE;
2435:     } else {
2436:       if (set) *set = PETSC_TRUE;
2437:       PetscOptionsStringToInt(value,ivalue);
2438:     }
2439:   } else {
2440:     if (set) *set = PETSC_FALSE;
2441:   }
2442:   return(0);
2443: }

2445: /*@C
2446:    PetscOptionsGetReal - Gets the double precision value for a particular
2447:    option in the database.

2449:    Not Collective

2451:    Input Parameters:
2452: +  options - options database, use NULL for default global database
2453: .  pre - string to prepend to each name or NULL
2454: -  name - the option one is seeking

2456:    Output Parameters:
2457: +  dvalue - the double value to return
2458: -  set - PETSC_TRUE if found, PETSC_FALSE if not found

2460:    Notes:
2461:     If the user does not supply the option dvalue is NOT changed. Thus
2462:      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2464:    Level: beginner

2466: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2467:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
2468:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2469:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2470:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2471:           PetscOptionsFList(), PetscOptionsEList()
2472: @*/
2473: PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set)
2474: {
2475:   const char     *value;
2476:   PetscBool      flag;

2482:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2483:   if (flag) {
2484:     if (!value) {
2485:       if (set) *set = PETSC_FALSE;
2486:     } else {
2487:       if (set) *set = PETSC_TRUE;
2488:       PetscOptionsStringToReal(value,dvalue);
2489:     }
2490:   } else {
2491:     if (set) *set = PETSC_FALSE;
2492:   }
2493:   return(0);
2494: }

2496: /*@C
2497:    PetscOptionsGetScalar - Gets the scalar value for a particular
2498:    option in the database.

2500:    Not Collective

2502:    Input Parameters:
2503: +  options - options database, use NULL for default global database
2504: .  pre - string to prepend to each name or NULL
2505: -  name - the option one is seeking

2507:    Output Parameters:
2508: +  dvalue - the double value to return
2509: -  set - PETSC_TRUE if found, else PETSC_FALSE

2511:    Level: beginner

2513:    Usage:
2514:    A complex number 2+3i must be specified with NO spaces

2516:    Notes:
2517:     If the user does not supply the option dvalue is NOT changed. Thus
2518:      you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true.

2520: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2521:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2522:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2523:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2524:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2525:           PetscOptionsFList(), PetscOptionsEList()
2526: @*/
2527: PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set)
2528: {
2529:   const char     *value;
2530:   PetscBool      flag;

2536:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2537:   if (flag) {
2538:     if (!value) {
2539:       if (set) *set = PETSC_FALSE;
2540:     } else {
2541: #if !defined(PETSC_USE_COMPLEX)
2542:       PetscOptionsStringToReal(value,dvalue);
2543: #else
2544:       PetscOptionsStringToScalar(value,dvalue);
2545: #endif
2546:       if (set) *set = PETSC_TRUE;
2547:     }
2548:   } else { /* flag */
2549:     if (set) *set = PETSC_FALSE;
2550:   }
2551:   return(0);
2552: }

2554: /*@C
2555:    PetscOptionsGetString - Gets the string value for a particular option in
2556:    the database.

2558:    Not Collective

2560:    Input Parameters:
2561: +  options - options database, use NULL for default global database
2562: .  pre - string to prepend to name or NULL
2563: .  name - the option one is seeking
2564: -  len - maximum length of the string including null termination

2566:    Output Parameters:
2567: +  string - location to copy string
2568: -  set - PETSC_TRUE if found, else PETSC_FALSE

2570:    Level: beginner

2572:    Fortran Note:
2573:    The Fortran interface is slightly different from the C/C++
2574:    interface (len is not used).  Sample usage in Fortran follows
2575: .vb
2576:       character *20    string
2577:       PetscErrorCode   ierr
2578:       PetscBool        set
2579:       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2580: .ve

2582:    Notes:
2583:     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

2585:            If the user does not use the option then the string is not changed. Thus
2586:            you should ALWAYS initialize the string if you access it without first checking if the set flag is true.

2588:     Note:
2589:       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).

2591: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2592:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2593:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2594:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2595:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2596:           PetscOptionsFList(), PetscOptionsEList()
2597: @*/
2598: PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set)
2599: {
2600:   const char     *value;
2601:   PetscBool      flag;

2607:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2608:   if (!flag) {
2609:     if (set) *set = PETSC_FALSE;
2610:   } else {
2611:     if (set) *set = PETSC_TRUE;
2612:     if (value) {
2613:       PetscStrncpy(string,value,len);
2614:     } else {
2615:       PetscArrayzero(string,len);
2616:     }
2617:   }
2618:   return(0);
2619: }

2621: char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[])
2622: {
2623:   const char     *value;
2624:   PetscBool      flag;

2628:   PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(NULL);
2629:   if (flag) PetscFunctionReturn((char*)value);
2630:   else PetscFunctionReturn(NULL);
2631: }

2633: /*@C
2634:    PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular
2635:    option in the database.  The values must be separated with commas with
2636:    no intervening spaces.

2638:    Not Collective

2640:    Input Parameters:
2641: +  options - options database, use NULL for default global database
2642: .  pre - string to prepend to each name or NULL
2643: -  name - the option one is seeking

2645:    Input/Output Parameter:
2646: .  nmax - maximum number of values to retrieve, on output the actual number of values retrieved

2648:    Output Parameters:
2649: +  dvalue - the integer values to return
2650: -  set - PETSC_TRUE if found, else PETSC_FALSE

2652:    Level: beginner

2654:    Notes:
2655:        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2656:        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE

2658: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2659:           PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2660:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2661:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2662:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2663:           PetscOptionsFList(), PetscOptionsEList()
2664: @*/
2665: PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set)
2666: {
2667:   const char     *svalue;
2668:   char           *value;
2670:   PetscInt       n = 0;
2671:   PetscBool      flag;
2672:   PetscToken     token;


2679:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2680:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2681:   if (set) *set = PETSC_TRUE;
2682:   PetscTokenCreate(svalue,',',&token);
2683:   PetscTokenFind(token,&value);
2684:   while (value && n < *nmax) {
2685:     PetscOptionsStringToBool(value,dvalue);
2686:     PetscTokenFind(token,&value);
2687:     dvalue++;
2688:     n++;
2689:   }
2690:   PetscTokenDestroy(&token);
2691:   *nmax = n;
2692:   return(0);
2693: }

2695: /*@C
2696:    PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database.

2698:    Not Collective

2700:    Input Parameters:
2701: +  options - options database, use NULL for default global database
2702: .  pre - option prefix or NULL
2703: .  name - option name
2704: -  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null

2706:    Input/Output Parameter:
2707: .  nmax - maximum number of values to retrieve, on output the actual number of values retrieved

2709:    Output Parameters:
2710: +  ivalue - the  enum values to return
2711: -  set - PETSC_TRUE if found, else PETSC_FALSE

2713:    Level: beginner

2715:    Notes:
2716:    The array must be passed as a comma separated list.

2718:    There must be no intervening spaces between the values.

2720:    list is usually something like PCASMTypes or some other predefined list of enum names.

2722: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2723:           PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2724:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(),
2725:           PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(),
2726:           PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2727:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2728: @*/
2729: PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set)
2730: {
2731:   const char     *svalue;
2732:   char           *value;
2733:   PetscInt       n = 0;
2734:   PetscEnum      evalue;
2735:   PetscBool      flag;
2736:   PetscToken     token;


2745:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2746:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2747:   if (set) *set = PETSC_TRUE;
2748:   PetscTokenCreate(svalue,',',&token);
2749:   PetscTokenFind(token,&value);
2750:   while (value && n < *nmax) {
2751:     PetscEnumFind(list,value,&evalue,&flag);
2752:     if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1);
2753:     ivalue[n++] = evalue;
2754:     PetscTokenFind(token,&value);
2755:   }
2756:   PetscTokenDestroy(&token);
2757:   *nmax = n;
2758:   return(0);
2759: }

2761: /*@C
2762:    PetscOptionsGetIntArray - Gets an array of integer values for a particular
2763:    option in the database.

2765:    Not Collective

2767:    Input Parameters:
2768: +  options - options database, use NULL for default global database
2769: .  pre - string to prepend to each name or NULL
2770: -  name - the option one is seeking

2772:    Input/Output Parameter:
2773: .  nmax - maximum number of values to retrieve, on output the actual number of values retrieved

2775:    Output Parameters:
2776: +  ivalue - the integer values to return
2777: -  set - PETSC_TRUE if found, else PETSC_FALSE

2779:    Level: beginner

2781:    Notes:
2782:    The array can be passed as
2783:    a comma separated list:                                 0,1,2,3,4,5,6,7
2784:    a range (start-end+1):                                  0-8
2785:    a range with given increment (start-end+1:inc):         0-7:2
2786:    a combination of values and ranges separated by commas: 0,1-8,8-15:2

2788:    There must be no intervening spaces between the values.

2790: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2791:           PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2792:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2793:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2794:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2795:           PetscOptionsFList(), PetscOptionsEList()
2796: @*/
2797: PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set)
2798: {
2799:   const char     *svalue;
2800:   char           *value;
2802:   PetscInt       n = 0,i,j,start,end,inc,nvalues;
2803:   size_t         len;
2804:   PetscBool      flag,foundrange;
2805:   PetscToken     token;


2812:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2813:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2814:   if (set) *set = PETSC_TRUE;
2815:   PetscTokenCreate(svalue,',',&token);
2816:   PetscTokenFind(token,&value);
2817:   while (value && n < *nmax) {
2818:     /* look for form  d-D where d and D are integers */
2819:     foundrange = PETSC_FALSE;
2820:     PetscStrlen(value,&len);
2821:     if (value[0] == '-') i=2;
2822:     else i=1;
2823:     for (;i<(int)len; i++) {
2824:       if (value[i] == '-') {
2825:         if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
2826:         value[i] = 0;

2828:         PetscOptionsStringToInt(value,&start);
2829:         inc  = 1;
2830:         j    = i+1;
2831:         for (;j<(int)len; j++) {
2832:           if (value[j] == ':') {
2833:             value[j] = 0;

2835:             PetscOptionsStringToInt(value+j+1,&inc);
2836:             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);
2837:             break;
2838:           }
2839:         }
2840:         PetscOptionsStringToInt(value+i+1,&end);
2841:         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);
2842:         nvalues = (end-start)/inc + (end-start)%inc;
2843:         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);
2844:         for (;start<end; start+=inc) {
2845:           *ivalue = start; ivalue++;n++;
2846:         }
2847:         foundrange = PETSC_TRUE;
2848:         break;
2849:       }
2850:     }
2851:     if (!foundrange) {
2852:       PetscOptionsStringToInt(value,ivalue);
2853:       ivalue++;
2854:       n++;
2855:     }
2856:     PetscTokenFind(token,&value);
2857:   }
2858:   PetscTokenDestroy(&token);
2859:   *nmax = n;
2860:   return(0);
2861: }

2863: /*@C
2864:    PetscOptionsGetRealArray - Gets an array of double precision values for a
2865:    particular option in the database.  The values must be separated with
2866:    commas with no intervening spaces.

2868:    Not Collective

2870:    Input Parameters:
2871: +  options - options database, use NULL for default global database
2872: .  pre - string to prepend to each name or NULL
2873: -  name - the option one is seeking

2875:    Input/Output Parameter:
2876: .  nmax - maximum number of values to retrieve, on output the actual number of values retrieved

2878:    Output Parameters:
2879: +  dvalue - the double values to return
2880: -  set - PETSC_TRUE if found, else PETSC_FALSE

2882:    Level: beginner

2884: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2885:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
2886:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2887:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2888:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2889:           PetscOptionsFList(), PetscOptionsEList()
2890: @*/
2891: PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set)
2892: {
2893:   const char     *svalue;
2894:   char           *value;
2896:   PetscInt       n = 0;
2897:   PetscBool      flag;
2898:   PetscToken     token;


2905:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2906:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2907:   if (set) *set = PETSC_TRUE;
2908:   PetscTokenCreate(svalue,',',&token);
2909:   PetscTokenFind(token,&value);
2910:   while (value && n < *nmax) {
2911:     PetscOptionsStringToReal(value,dvalue++);
2912:     PetscTokenFind(token,&value);
2913:     n++;
2914:   }
2915:   PetscTokenDestroy(&token);
2916:   *nmax = n;
2917:   return(0);
2918: }

2920: /*@C
2921:    PetscOptionsGetScalarArray - Gets an array of scalars for a
2922:    particular option in the database.  The values must be separated with
2923:    commas with no intervening spaces.

2925:    Not Collective

2927:    Input Parameters:
2928: +  options - options database, use NULL for default global database
2929: .  pre - string to prepend to each name or NULL
2930: -  name - the option one is seeking

2932:    Input/Output Parameter:
2933: .  nmax - maximum number of values to retrieve, on output the actual number of values retrieved

2935:    Output Parameters:
2936: +  dvalue - the scalar values to return
2937: -  set - PETSC_TRUE if found, else PETSC_FALSE

2939:    Level: beginner

2941: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2942:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
2943:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2944:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2945:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2946:           PetscOptionsFList(), PetscOptionsEList()
2947: @*/
2948: PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],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:     PetscOptionsStringToScalar(value,dvalue++);
2969:     PetscTokenFind(token,&value);
2970:     n++;
2971:   }
2972:   PetscTokenDestroy(&token);
2973:   *nmax = n;
2974:   return(0);
2975: }

2977: /*@C
2978:    PetscOptionsGetStringArray - Gets an array of string values for a particular
2979:    option in the database. The values must be separated with commas with
2980:    no intervening spaces.

2982:    Not Collective

2984:    Input Parameters:
2985: +  options - options database, use NULL for default global database
2986: .  pre - string to prepend to name or NULL
2987: -  name - the option one is seeking

2989:    Input/Output Parameter:
2990: .  nmax - maximum number of strings, on output the actual number of strings found

2992:    Output Parameters:
2993: +  strings - location to copy strings
2994: -  set - PETSC_TRUE if found, else PETSC_FALSE

2996:    Level: beginner

2998:    Notes:
2999:    The nmax parameter is used for both input and output.

3001:    The user should pass in an array of pointers to char, to hold all the
3002:    strings returned by this function.

3004:    The user is responsible for deallocating the strings that are
3005:    returned. The Fortran interface for this routine is not supported.

3007: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
3008:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
3009:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
3010:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
3011:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
3012:           PetscOptionsFList(), PetscOptionsEList()
3013: @*/
3014: PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set)
3015: {
3016:   const char     *svalue;
3017:   char           *value;
3019:   PetscInt       n = 0;
3020:   PetscBool      flag;
3021:   PetscToken     token;


3028:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
3029:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
3030:   if (set) *set = PETSC_TRUE;
3031:   PetscTokenCreate(svalue,',',&token);
3032:   PetscTokenFind(token,&value);
3033:   while (value && n < *nmax) {
3034:     PetscStrallocpy(value,&strings[n]);
3035:     PetscTokenFind(token,&value);
3036:     n++;
3037:   }
3038:   PetscTokenDestroy(&token);
3039:   *nmax = n;
3040:   return(0);
3041: }

3043: /*@C
3044:    PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one

3046:    Prints a deprecation warning, unless an option is supplied to suppress.

3048:    Logically Collective

3050:    Input Parameters:
3051: +  pre - string to prepend to name or NULL
3052: .  oldname - the old, deprecated option
3053: .  newname - the new option, or NULL if option is purely removed
3054: .  version - a string describing the version of first deprecation, e.g. "3.9"
3055: -  info - additional information string, or NULL.

3057:    Options Database Keys:
3058: . -options_suppress_deprecated_warnings - do not print deprecation warnings

3060:    Notes:
3061:    Must be called between PetscOptionsBegin() (or PetscObjectOptionsBegin()) and PetscOptionsEnd().
3062:    Only the proces of rank zero that owns the PetscOptionsItems are argument (managed by PetscOptionsBegin() or
3063:    PetscObjectOptionsBegin() prints the information
3064:    If newname is provided, the old option is replaced. Otherwise, it remains
3065:    in the options database.
3066:    If an option is not replaced, the info argument should be used to advise the user
3067:    on how to proceed.
3068:    There is a limit on the length of the warning printed, so very long strings
3069:    provided as info may be truncated.

3071:    Level: developer

3073: .seealso: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsScalar(), PetscOptionsBool(), PetscOptionsString(), PetscOptionsSetValue()

3075: @*/
3076: PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[])
3077: {
3078:   PetscErrorCode     ierr;
3079:   PetscBool          found,quiet;
3080:   const char         *value;
3081:   const char * const quietopt="-options_suppress_deprecated_warnings";
3082:   char               msg[4096];
3083:   char               *prefix = NULL;
3084:   PetscOptions       options = NULL;
3085:   MPI_Comm           comm = PETSC_COMM_SELF;

3090:   if (PetscOptionsObject) {
3091:     prefix  = PetscOptionsObject->prefix;
3092:     options = PetscOptionsObject->options;
3093:     comm    = PetscOptionsObject->comm;
3094:   }
3095:   PetscOptionsFindPair(options,prefix,oldname,&value,&found);
3096:   if (found) {
3097:     if (newname) {
3098:       if (prefix) {
3099:         PetscOptionsPrefixPush(options,prefix);
3100:       }
3101:       PetscOptionsSetValue(options,newname,value);
3102:       if (prefix) {
3103:         PetscOptionsPrefixPop(options);
3104:       }
3105:       PetscOptionsClearValue(options,oldname);
3106:     }
3107:     quiet = PETSC_FALSE;
3108:     PetscOptionsGetBool(options,NULL,quietopt,&quiet,NULL);
3109:     if (!quiet) {
3110:       PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");
3111:       PetscStrcat(msg,oldname);
3112:       PetscStrcat(msg," is deprecated as of version ");
3113:       PetscStrcat(msg,version);
3114:       PetscStrcat(msg," and will be removed in a future release.");
3115:       if (newname) {
3116:         PetscStrcat(msg," Please use the option ");
3117:         PetscStrcat(msg,newname);
3118:         PetscStrcat(msg," instead.");
3119:       }
3120:       if (info) {
3121:         PetscStrcat(msg," ");
3122:         PetscStrcat(msg,info);
3123:       }
3124:       PetscStrcat(msg," (Silence this warning with ");
3125:       PetscStrcat(msg,quietopt);
3126:       PetscStrcat(msg,")\n");
3127:       PetscPrintf(comm,msg);
3128:     }
3129:   }
3130:   return(0);
3131: }