Actual source code: options.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: /* Define Feature test macros to make sure atoll is available (SVr4, POSIX.1-2001, 4.3BSD, C99), not in (C89 and POSIX.1-1996) */
  2: #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for atoll() */

  4: /*
  5:    These routines simplify the use of command line, file options, etc., and are used to manipulate the options database.
  6:    This provides the low-level interface, the high level interface is in aoptions.c

  8:    Some routines use regular malloc and free because it cannot know  what malloc is requested with the
  9:    options database until it has already processed the input.
 10: */

 12: #include <petsc/private/petscimpl.h>
 13: #include <petscviewer.h>
 14: #include <ctype.h>
 15: #if defined(PETSC_HAVE_MALLOC_H)
 16: #include <malloc.h>
 17: #endif
 18: #if defined(PETSC_HAVE_STRINGS_H)
 19: #  include <strings.h>          /* strcasecmp */
 20: #endif
 21: #if defined(PETSC_HAVE_YAML)
 22: #include <yaml.h>
 23: #endif

 25: #if defined(PETSC_HAVE_STRCASECMP)
 26: #define PetscOptNameCmp(a,b) strcasecmp(a,b)
 27: #elif defined(PETSC_HAVE_STRICMP)
 28: #define PetscOptNameCmp(a,b) stricmp(a,b)
 29: #else
 30: #define PetscOptNameCmp(a,b) Error_strcasecmp_not_found
 31: #endif

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

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

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

 57: PETSC_STATIC_INLINE int PetscOptEqual(const char a[],const char b[])
 58: {
 59:   return !PetscOptNameCmp(a,b);
 60: }

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

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

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

 81:   /* Hash table */
 82:   khash_t(HO)    *ht;

 84:   /* Prefixes */
 85:   int            prefixind;
 86:   int            prefixstack[MAXPREFIXES];
 87:   char           prefix[MAXOPTNAME];

 89:   /* Aliases */
 90:   int            Naliases;                   /* number or aliases */
 91:   char           *aliases1[MAXALIASES];      /* aliased */
 92:   char           *aliases2[MAXALIASES];      /* aliasee */

 94:   /* Help */
 95:   PetscBool      help;       /* flag whether "-help" is in the database */
 96:   PetscBool      help_intro; /* flag whether "-help intro" is in the database */

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

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

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

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

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

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

134: /*@
135:    PetscOptionsCreate - Creates an empty options database.

137:    Logically collective

139:    Output Parameter:
140: .  options - Options database object

142:    Level: advanced

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

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

156: /*@
157:     PetscOptionsDestroy - Destroys an option database.

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

161:   Input Parameter:
162: .  options - the PetscOptions object

164:    Level: advanced

166: .seealso: PetscOptionsInsert(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsInsert(), PetscOptionsSetValue()
167: @*/
168: PetscErrorCode PetscOptionsDestroy(PetscOptions *options)
169: {

172:   if (!*options) return 0;
173:   if ((*options)->previous) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You are destroying an option that has been used with PetscOptionsPush() but does not have a corresponding PetscOptionsPop()");
174:   PetscOptionsClear(*options);if (ierr) return ierr;
175:   /* XXX what about monitors ? */
176:   free(*options);
177:   *options = NULL;
178:   return(0);
179: }

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

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

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

198:   Logically Collective

200:   Input Parameter:
201: .   opt - the options obtained with PetscOptionsCreate()

203:   Notes:
204:   Use PetscOptionsPop() to return to the previous default options database

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

211:    Level: advanced

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

215: @*/
216: PetscErrorCode PetscOptionsPush(PetscOptions opt)
217: {

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

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

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

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

236:    Level: advanced

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

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

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

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

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

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

275:    Not collective

277:    Input Parameter:
278: .  key - string to check if valid

280:    Output Parameter:
281: .  valid - PETSC_TRUE if a valid key

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

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

303: /*@C
304:    PetscOptionsInsertString - Inserts options into the database from a string

306:    Logically Collective

308:    Input Parameter:
309: +  options - options object
310: -  in_str - string that contains options separated by blanks

312:    Level: intermediate

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

319:    Contributed by Boyana Norris

321: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
322:           PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
323:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
324:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
325:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
326:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsInsertFile()
327: @*/
328: PetscErrorCode PetscOptionsInsertString(PetscOptions options,const char in_str[])
329: {
330:   char           *first,*second;
332:   PetscToken     token;
333:   PetscBool      key,ispush,ispop,isopts;

336:   PetscTokenCreate(in_str,' ',&token);
337:   PetscTokenFind(token,&first);
338:   while (first) {
339:     PetscStrcasecmp(first,"-prefix_push",&ispush);
340:     PetscStrcasecmp(first,"-prefix_pop",&ispop);
341:     PetscStrcasecmp(first,"-options_file",&isopts);
342:     PetscOptionsValidKey(first,&key);
343:     if (ispush) {
344:       PetscTokenFind(token,&second);
345:       PetscOptionsPrefixPush(options,second);
346:       PetscTokenFind(token,&first);
347:     } else if (ispop) {
348:       PetscOptionsPrefixPop(options);
349:       PetscTokenFind(token,&first);
350:     } else if (isopts) {
351:       PetscTokenFind(token,&second);
352:       PetscOptionsInsertFile(PETSC_COMM_SELF,options,second,PETSC_TRUE);
353:       PetscTokenFind(token,&first);
354:     } else if (key) {
355:       PetscTokenFind(token,&second);
356:       PetscOptionsValidKey(second,&key);
357:       if (!key) {
358:         PetscOptionsSetValue(options,first,second);
359:         PetscTokenFind(token,&first);
360:       } else {
361:         PetscOptionsSetValue(options,first,NULL);
362:         first = second;
363:       }
364:     } else {
365:       PetscTokenFind(token,&first);
366:     }
367:   }
368:   PetscTokenDestroy(&token);
369:   return(0);
370: }

372: /*
373:     Returns a line (ended by a \n, \r or null character of any length. Result should be freed with free()
374: */
375: static char *Petscgetline(FILE * f)
376: {
377:   size_t size  = 0;
378:   size_t len   = 0;
379:   size_t last  = 0;
380:   char   *buf  = NULL;

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

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

400:      Collective

402:   Input Parameter:
403: +   comm - the processes that will share the options (usually PETSC_COMM_WORLD)
404: .   options - options database, use NULL for default global database
405: .   file - name of file
406: -   require - if PETSC_TRUE will generate an error if the file does not exist


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

419:   Level: developer

421: .seealso: PetscOptionsSetValue(), PetscOptionsView(), PetscOptionsHasName(), PetscOptionsGetInt(),
422:           PetscOptionsGetReal(), PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
423:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
424:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
425:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
426:           PetscOptionsFList(), PetscOptionsEList()

428: @*/
429: PetscErrorCode PetscOptionsInsertFile(MPI_Comm comm,PetscOptions options,const char file[],PetscBool require)
430: {
431:   char           *string,fname[PETSC_MAX_PATH_LEN],*vstring = NULL,*astring = NULL,*packed = NULL;
432:   char           *tokens[4];
434:   size_t         i,len,bytes;
435:   FILE           *fd;
436:   PetscToken     token=NULL;
437:   int            err;
438:   char           *cmatch;
439:   const char     cmt='#';
440:   PetscInt       line=1;
441:   PetscMPIInt    rank,cnt=0,acnt=0,counts[2];
442:   PetscBool      isdir,alias=PETSC_FALSE,valid;

445:   MPI_Comm_rank(comm,&rank);
446:   PetscMemzero(tokens,sizeof(tokens));
447:   if (!rank) {
448:     cnt        = 0;
449:     acnt       = 0;

451:     PetscFixFilename(file,fname);
452:     fd   = fopen(fname,"r");
453:     PetscTestDirectory(fname,'r',&isdir);
454:     if (isdir && require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified options file %s is a directory",fname);
455:     if (fd && !isdir) {
456:       PetscSegBuffer vseg,aseg;
457:       PetscSegBufferCreate(1,4000,&vseg);
458:       PetscSegBufferCreate(1,2000,&aseg);

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

463:       while ((string = Petscgetline(fd))) {
464:         /* eliminate comments from each line */
465:         PetscStrchr(string,cmt,&cmatch);
466:         if (cmatch) *cmatch = 0;
467:         PetscStrlen(string,&len);
468:         /* replace tabs, ^M, \n with " " */
469:         for (i=0; i<len; i++) {
470:           if (string[i] == '\t' || string[i] == '\r' || string[i] == '\n') {
471:             string[i] = ' ';
472:           }
473:         }
474:         PetscTokenCreate(string,' ',&token);
475:         PetscTokenFind(token,&tokens[0]);
476:         if (!tokens[0]) {
477:           goto destroy;
478:         } else if (!tokens[0][0]) { /* if token 0 is empty (string begins with spaces), redo */
479:           PetscTokenFind(token,&tokens[0]);
480:         }
481:         for (i=1; i<4; i++) {
482:           PetscTokenFind(token,&tokens[i]);
483:         }
484:         if (!tokens[0]) {
485:           goto destroy;
486:         } else if (tokens[0][0] == '-') {
487:           PetscOptionsValidKey(tokens[0],&valid);
488:           if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: invalid option %s",fname,line,tokens[0]);
489:           PetscStrlen(tokens[0],&len);
490:           PetscSegBufferGet(vseg,len+1,&vstring);
491:           PetscArraycpy(vstring,tokens[0],len);
492:           vstring[len] = ' ';
493:           if (tokens[1]) {
494:             PetscOptionsValidKey(tokens[1],&valid);
495:             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]);
496:             PetscStrlen(tokens[1],&len);
497:             PetscSegBufferGet(vseg,len+3,&vstring);
498:             vstring[0] = '"';
499:             PetscArraycpy(vstring+1,tokens[1],len);
500:             vstring[len+1] = '"';
501:             vstring[len+2] = ' ';
502:           }
503:         } else {
504:           PetscStrcasecmp(tokens[0],"alias",&alias);
505:           if (alias) {
506:             PetscOptionsValidKey(tokens[1],&valid);
507:             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]);
508:             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]);
509:             PetscOptionsValidKey(tokens[2],&valid);
510:             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]);
511:             PetscStrlen(tokens[1],&len);
512:             PetscSegBufferGet(aseg,len+1,&astring);
513:             PetscArraycpy(astring,tokens[1],len);
514:             astring[len] = ' ';

516:             PetscStrlen(tokens[2],&len);
517:             PetscSegBufferGet(aseg,len+1,&astring);
518:             PetscArraycpy(astring,tokens[2],len);
519:             astring[len] = ' ';
520:           } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown first token in options file %s line %D: %s",fname,line,tokens[0]);
521:         }
522:         {
523:           const char *extraToken = alias ? tokens[3] : tokens[2];
524:           if (extraToken) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in options file %s line %D: extra token %s",fname,line,extraToken);
525:         }
526: destroy:
527:         free(string);
528:         PetscTokenDestroy(&token);
529:         alias = PETSC_FALSE;
530:         line++;
531:       }
532:       err = fclose(fd);
533:       if (err) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file %s",fname);
534:       PetscSegBufferGetSize(aseg,&bytes); /* size without null termination */
535:       PetscMPIIntCast(bytes,&acnt);
536:       PetscSegBufferGet(aseg,1,&astring);
537:       astring[0] = 0;
538:       PetscSegBufferGetSize(vseg,&bytes); /* size without null termination */
539:       PetscMPIIntCast(bytes,&cnt);
540:       PetscSegBufferGet(vseg,1,&vstring);
541:       vstring[0] = 0;
542:       PetscMalloc1(2+acnt+cnt,&packed);
543:       PetscSegBufferExtractTo(aseg,packed);
544:       PetscSegBufferExtractTo(vseg,packed+acnt+1);
545:       PetscSegBufferDestroy(&aseg);
546:       PetscSegBufferDestroy(&vseg);
547:     } else if (require) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Unable to open options file %s",fname);
548:   }

550:   counts[0] = acnt;
551:   counts[1] = cnt;
552:   err = MPI_Bcast(counts,2,MPI_INT,0,comm);
553:   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://www.mcs.anl.gov/petsc/documentation/faq.html");
554:   acnt = counts[0];
555:   cnt = counts[1];
556:   if (rank) {
557:     PetscMalloc1(2+acnt+cnt,&packed);
558:   }
559:   if (acnt || cnt) {
560:     MPI_Bcast(packed,2+acnt+cnt,MPI_CHAR,0,comm);
561:     astring = packed;
562:     vstring = packed + acnt + 1;
563:   }

565:   if (acnt) {
566:     PetscTokenCreate(astring,' ',&token);
567:     PetscTokenFind(token,&tokens[0]);
568:     while (tokens[0]) {
569:       PetscTokenFind(token,&tokens[1]);
570:       PetscOptionsSetAlias(options,tokens[0],tokens[1]);
571:       PetscTokenFind(token,&tokens[0]);
572:     }
573:     PetscTokenDestroy(&token);
574:   }

576:   if (cnt) {
577:     PetscOptionsInsertString(options,vstring);
578:   }
579:   PetscFree(packed);
580:   return(0);
581: }

583: static PetscErrorCode PetscOptionsInsertArgs(PetscOptions options,int argc,char *args[])
584: {
586:   int            left    = argc - 1;
587:   char           **eargs = args + 1;

590:   while (left) {
591:     PetscBool isoptions_file,isprefixpush,isprefixpop,isp4,tisp4,isp4yourname,isp4rmrank,key;
592:     PetscStrcasecmp(eargs[0],"-options_file",&isoptions_file);
593:     PetscStrcasecmp(eargs[0],"-prefix_push",&isprefixpush);
594:     PetscStrcasecmp(eargs[0],"-prefix_pop",&isprefixpop);
595:     PetscStrcasecmp(eargs[0],"-p4pg",&isp4);
596:     PetscStrcasecmp(eargs[0],"-p4yourname",&isp4yourname);
597:     PetscStrcasecmp(eargs[0],"-p4rmrank",&isp4rmrank);
598:     PetscStrcasecmp(eargs[0],"-p4wd",&tisp4);
599:     isp4 = (PetscBool) (isp4 || tisp4);
600:     PetscStrcasecmp(eargs[0],"-np",&tisp4);
601:     isp4 = (PetscBool) (isp4 || tisp4);
602:     PetscStrcasecmp(eargs[0],"-p4amslave",&tisp4);
603:     PetscOptionsValidKey(eargs[0],&key);

605:     if (!key) {
606:       eargs++; left--;
607:     } else if (isoptions_file) {
608:       if (left <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file filename option");
609:       if (eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing filename for -options_file filename option");
610:       PetscOptionsInsertFile(PETSC_COMM_WORLD,options,eargs[1],PETSC_TRUE);
611:       eargs += 2; left -= 2;
612:     } else if (isprefixpush) {
613:       if (left <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option");
614:       if (eargs[1][0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing prefix for -prefix_push option (prefixes cannot start with '-')");
615:       PetscOptionsPrefixPush(options,eargs[1]);
616:       eargs += 2; left -= 2;
617:     } else if (isprefixpop) {
618:       PetscOptionsPrefixPop(options);
619:       eargs++; left--;

621:       /*
622:        These are "bad" options that MPICH, etc put on the command line
623:        we strip them out here.
624:        */
625:     } else if (tisp4 || isp4rmrank) {
626:       eargs += 1; left -= 1;
627:     } else if (isp4 || isp4yourname) {
628:       eargs += 2; left -= 2;
629:     } else {
630:       PetscBool nextiskey = PETSC_FALSE;
631:       if (left >= 2) {PetscOptionsValidKey(eargs[1],&nextiskey);}
632:       if (left < 2 || nextiskey) {
633:         PetscOptionsSetValue(options,eargs[0],NULL);
634:         eargs++; left--;
635:       } else {
636:         PetscOptionsSetValue(options,eargs[0],eargs[1]);
637:         eargs += 2; left -= 2;
638:       }
639:     }
640:   }
641:   return(0);
642: }

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

649:   if (set[opt]) {
650:     PetscOptionsStringToBool(val[opt],flg);
651:   } else *flg = PETSC_FALSE;
652:   return(0);
653: }

655: /* Process options with absolute precedence */
656: static PetscErrorCode PetscOptionsProcessPrecedentFlags(PetscOptions options,int argc,char *args[],PetscBool *skip_petscrc,PetscBool *skip_petscrc_set)
657: {
658:   const char* const *opt = precedentOptions;
659:   const size_t      n = PO_NUM;
660:   size_t            o;
661:   int               a;
662:   const char        **val;
663:   PetscBool         *set;
664:   PetscErrorCode    ierr;

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

669:   /* Look for options possibly set using PetscOptionsSetValue beforehand */
670:   for (o=0; o<n; o++) {
671:     PetscOptionsFindPair(options,NULL,opt[o],&val[o],&set[o]);
672:   }

674:   /* Loop through all args to collect last occuring value of each option */
675:   for (a=1; a<argc; a++) {
676:     PetscBool valid, eq;

678:     PetscOptionsValidKey(args[a],&valid);
679:     if (!valid) continue;
680:     for (o=0; o<n; o++) {
681:       PetscStrcasecmp(args[a],opt[o],&eq);
682:       if (eq) {
683:         set[o] = PETSC_TRUE;
684:         if (a == argc-1 || !args[a+1] || !args[a+1][0] || args[a+1][0] == '-') val[o] = NULL;
685:         else val[o] = args[a+1];
686:         break;
687:       }
688:     }
689:   }

691:   /* Process flags */
692:   PetscStrcasecmp(val[PO_HELP], "intro", &options->help_intro);
693:   if (options->help_intro) options->help = PETSC_TRUE;
694:   else {PetscOptionsStringToBoolIfSet_Private(PO_HELP,            val,set,&options->help);}
695:   PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR_CANCEL,val,set,&options->monitorCancel);
696:   PetscOptionsStringToBoolIfSet_Private(PO_OPTIONS_MONITOR,       val,set,&options->monitorFromOptions);
697:   PetscOptionsStringToBoolIfSet_Private(PO_SKIP_PETSCRC,          val,set,skip_petscrc);
698:   *skip_petscrc_set = set[PO_SKIP_PETSCRC];

700:   /* Store precedent options in database and mark them as used */
701:   for (o=0; o<n; o++) {
702:     if (set[o]) {
703:       int pos;

705:       PetscOptionsSetValue_Private(options,opt[o],val[o],&pos);
706:       options->used[pos] = PETSC_TRUE;
707:     }
708:   }

710:   PetscFree2(val,set);
711:   options->precedentProcessed = PETSC_TRUE;
712:   return(0);
713: }

715: PETSC_STATIC_INLINE PetscErrorCode PetscOptionsSkipPrecedent(PetscOptions options,const char name[],PetscBool *flg)
716: {
717:   int i;

720:   *flg = PETSC_FALSE;
721:   if (options->precedentProcessed) {
722:     for (i=0; i<PO_NUM; i++) {
723:       if (!PetscOptNameCmp(precedentOptions[i],name)) {
724:         /* check if precedent option has been set already */
725:         PetscOptionsFindPair(options,NULL,name,NULL,flg);
726:         if (*flg) break;
727:       }
728:     }
729:   }
730:   return(0);
731: }

733: /*@C
734:    PetscOptionsInsert - Inserts into the options database from the command line,
735:                         the environmental variable and a file.

737:    Collective on PETSC_COMM_WORLD

739:    Input Parameters:
740: +  options - options database or NULL for the default global database
741: .  argc - count of number of command line arguments
742: .  args - the command line arguments
743: -  file - [optional] PETSc database file, also checks ~/.petscrc, .petscrc and petscrc.
744:           Use NULL to not check for code specific file.
745:           Use -skip_petscrc in the code specific file (or command line) to skip ~/.petscrc, .petscrc and petscrc files.

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

752:    Options Database Keys:
753: .   -options_file <filename> - read options from a file

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

757:    Level: advanced

759: .seealso: PetscOptionsDestroy(), PetscOptionsView(), PetscOptionsInsertString(), PetscOptionsInsertFile(),
760:           PetscInitialize()
761: @*/
762: PetscErrorCode PetscOptionsInsert(PetscOptions options,int *argc,char ***args,const char file[])
763: {
765:   PetscMPIInt    rank;
766:   char           filename[PETSC_MAX_PATH_LEN];
767:   PetscBool      hasArgs = (argc && *argc) ? PETSC_TRUE : PETSC_FALSE;
768:   PetscBool      skipPetscrc = PETSC_FALSE, skipPetscrcSet = PETSC_FALSE;


772:   if (hasArgs && !(args && *args)) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_NULL, "*argc > 1 but *args not given");
773:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

775:   if (!options) {
776:     PetscOptionsCreateDefault();
777:     options = defaultoptions;
778:   }
779:   if (hasArgs) {
780:     /* process options with absolute precedence */
781:     PetscOptionsProcessPrecedentFlags(options,*argc,*args,&skipPetscrc,&skipPetscrcSet);
782:   }
783:   if (file && file[0]) {
784:     PetscStrreplace(PETSC_COMM_WORLD,file,filename,PETSC_MAX_PATH_LEN);
785:     PetscOptionsInsertFile(PETSC_COMM_WORLD,options,filename,PETSC_TRUE);
786:     /* if -skip_petscrc has not been set from command line, check whether it has been set in the file */
787:     if (!skipPetscrcSet) {PetscOptionsGetBool(options,NULL,"-skip_petscrc",&skipPetscrc,NULL);}
788:   }
789:   if (!skipPetscrc) {
790:     PetscGetHomeDirectory(filename,PETSC_MAX_PATH_LEN-16);
791:     /* PetscOptionsInsertFile() does a fopen() on rank0 only - so only rank0 HomeDir value is relavent */
792:     if (filename[0]) { PetscStrcat(filename,"/.petscrc"); }
793:     PetscOptionsInsertFile(PETSC_COMM_WORLD,options,filename,PETSC_FALSE);
794:     PetscOptionsInsertFile(PETSC_COMM_WORLD,options,".petscrc",PETSC_FALSE);
795:     PetscOptionsInsertFile(PETSC_COMM_WORLD,options,"petscrc",PETSC_FALSE);
796:   }

798:   /* insert environment options */
799:   {
800:     char   *eoptions = NULL;
801:     size_t len       = 0;
802:     if (!rank) {
803:       eoptions = (char*)getenv("PETSC_OPTIONS");
804:       PetscStrlen(eoptions,&len);
805:       MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
806:     } else {
807:       MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
808:       if (len) {
809:         PetscMalloc1(len+1,&eoptions);
810:       }
811:     }
812:     if (len) {
813:       MPI_Bcast(eoptions,len,MPI_CHAR,0,PETSC_COMM_WORLD);
814:       if (rank) eoptions[len] = 0;
815:       PetscOptionsInsertString(options,eoptions);
816:       if (rank) {PetscFree(eoptions);}
817:     }
818:   }

820: #if defined(PETSC_HAVE_YAML)
821:   {
822:     char   *eoptions = NULL;
823:     size_t len       = 0;
824:     if (!rank) {
825:       eoptions = (char*)getenv("PETSC_OPTIONS_YAML");
826:       PetscStrlen(eoptions,&len);
827:       MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
828:     } else {
829:       MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);
830:       if (len) {
831:         PetscMalloc1(len+1,&eoptions);
832:       }
833:     }
834:     if (len) {
835:       MPI_Bcast(eoptions,len,MPI_CHAR,0,PETSC_COMM_WORLD);
836:       if (rank) eoptions[len] = 0;
837:       PetscOptionsInsertStringYAML(options,eoptions);
838:       if (rank) {PetscFree(eoptions);}
839:     }
840:   }
841:   {
842:     char      yaml_file[PETSC_MAX_PATH_LEN];
843:     char      yaml_string[BUFSIZ];
844:     PetscBool yaml_flg;
845:     PetscOptionsGetString(NULL,NULL,"-options_file_yaml",yaml_file,sizeof(yaml_file),&yaml_flg);
846:     if (yaml_flg) {
847:       PetscOptionsInsertFileYAML(PETSC_COMM_WORLD,yaml_file,PETSC_TRUE);
848:     }
849:     PetscOptionsGetString(NULL,NULL,"-options_string_yaml",yaml_string,sizeof(yaml_string),&yaml_flg);
850:     if (yaml_flg) {
851:       PetscOptionsInsertStringYAML(NULL,yaml_string);
852:     }
853:   }
854: #endif

856:   /* insert command line options here because they take precedence over arguments in petscrc/environment */
857:   if (hasArgs) {PetscOptionsInsertArgs(options,*argc,*args);}
858:   return(0);
859: }

861: /*@C
862:    PetscOptionsView - Prints the options that have been loaded. This is
863:    useful for debugging purposes.

865:    Logically Collective on PetscViewer

867:    Input Parameter:
868: +  options - options database, use NULL for default global database
869: -  viewer - must be an PETSCVIEWERASCII viewer

871:    Options Database Key:
872: .  -options_view - Activates PetscOptionsView() within PetscFinalize()

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

878:    Level: advanced

880: .seealso: PetscOptionsAllUsed()
881: @*/
882: PetscErrorCode PetscOptionsView(PetscOptions options,PetscViewer viewer)
883: {
885:   PetscInt       i;
886:   PetscBool      isascii;

890:   options = options ? options : defaultoptions;
891:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD;
892:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
893:   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Only supports ASCII viewer");

895:   if (!options->N) {
896:     PetscViewerASCIIPrintf(viewer,"#No PETSc Option Table entries\n");
897:     return(0);
898:   }

900:   PetscViewerASCIIPrintf(viewer,"#PETSc Option Table entries:\n");
901:   for (i=0; i<options->N; i++) {
902:     if (options->values[i]) {
903:       PetscViewerASCIIPrintf(viewer,"-%s %s\n",options->names[i],options->values[i]);
904:     } else {
905:       PetscViewerASCIIPrintf(viewer,"-%s\n",options->names[i]);
906:     }
907:   }
908:   PetscViewerASCIIPrintf(viewer,"#End of PETSc Option Table entries\n");
909:   return(0);
910: }

912: /*
913:    Called by error handlers to print options used in run
914: */
915: PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void)
916: {
917:   PetscInt     i;
918:   PetscOptions options = defaultoptions;

921:   if (options->N) {
922:     (*PetscErrorPrintf)("PETSc Option Table entries:\n");
923:   } else {
924:     (*PetscErrorPrintf)("No PETSc Option Table entries\n");
925:   }
926:   for (i=0; i<options->N; i++) {
927:     if (options->values[i]) {
928:       (*PetscErrorPrintf)("-%s %s\n",options->names[i],options->values[i]);
929:     } else {
930:       (*PetscErrorPrintf)("-%s\n",options->names[i]);
931:     }
932:   }
933:   return(0);
934: }

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

939:    Logically Collective

941:    Input Parameter:
942: +  options - options database, or NULL for the default global database
943: -  prefix - The string to append to the existing prefix

945:    Options Database Keys:
946: +   -prefix_push <some_prefix_> - push the given prefix
947: -   -prefix_pop - pop the last prefix

949:    Notes:
950:    It is common to use this in conjunction with -options_file as in

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

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

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

961: Level: advanced

963: .seealso: PetscOptionsPrefixPop(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
964: @*/
965: PetscErrorCode PetscOptionsPrefixPush(PetscOptions options,const char prefix[])
966: {
968:   size_t         n;
969:   PetscInt       start;
970:   char           key[MAXOPTNAME+1];
971:   PetscBool      valid;

975:   options = options ? options : defaultoptions;
976:   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);
977:   key[0] = '-'; /* keys must start with '-' */
978:   PetscStrncpy(key+1,prefix,sizeof(key)-1);
979:   PetscOptionsValidKey(key,&valid);
980:   if (!valid && options->prefixind > 0 && isdigit((int)prefix[0])) valid = PETSC_TRUE; /* If the prefix stack is not empty, make numbers a valid prefix */
981:   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":"");
982:   start = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
983:   PetscStrlen(prefix,&n);
984:   if (n+1 > sizeof(options->prefix)-start) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Maximum prefix length %d exceeded",sizeof(options->prefix));
985:   PetscArraycpy(options->prefix+start,prefix,n+1);
986:   options->prefixstack[options->prefixind++] = start+n;
987:   return(0);
988: }

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

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

995:   Input Parameters:
996: .  options - options database, or NULL for the default global database

998:    Level: advanced

1000: .seealso: PetscOptionsPrefixPush(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue()
1001: @*/
1002: PetscErrorCode PetscOptionsPrefixPop(PetscOptions options)
1003: {
1004:   PetscInt offset;

1007:   options = options ? options : defaultoptions;
1008:   if (options->prefixind < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"More prefixes popped than pushed");
1009:   options->prefixind--;
1010:   offset = options->prefixind ? options->prefixstack[options->prefixind-1] : 0;
1011:   options->prefix[offset] = 0;
1012:   return(0);
1013: }

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

1018:     Logically Collective

1020:   Input Parameters:
1021: .  options - options database, use NULL for the default global database

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

1028:    Level: developer

1030: .seealso: PetscOptionsInsert()
1031: @*/
1032: PetscErrorCode PetscOptionsClear(PetscOptions options)
1033: {
1034:   PetscInt i;

1036:   options = options ? options : defaultoptions;
1037:   if (!options) return 0;

1039:   for (i=0; i<options->N; i++) {
1040:     if (options->names[i])  free(options->names[i]);
1041:     if (options->values[i]) free(options->values[i]);
1042:   }
1043:   options->N = 0;

1045:   for (i=0; i<options->Naliases; i++) {
1046:     free(options->aliases1[i]);
1047:     free(options->aliases2[i]);
1048:   }
1049:   options->Naliases = 0;

1051:   /* destroy hash table */
1052:   kh_destroy(HO,options->ht);
1053:   options->ht = NULL;

1055:   options->prefixind = 0;
1056:   options->prefix[0] = 0;
1057:   options->help      = PETSC_FALSE;
1058:   return 0;
1059: }

1061: /*@C
1062:    PetscOptionsSetAlias - Makes a key and alias for another key

1064:    Logically Collective

1066:    Input Parameters:
1067: +  options - options database, or NULL for default global database
1068: .  newname - the alias
1069: -  oldname - the name that alias will refer to

1071:    Level: advanced

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

1078: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1079:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
1080:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1081:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1082:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1083:           PetscOptionsFList(), PetscOptionsEList()
1084: @*/
1085: PetscErrorCode PetscOptionsSetAlias(PetscOptions options,const char newname[],const char oldname[])
1086: {
1087:   PetscInt       n;
1088:   size_t         len;
1089:   PetscBool      valid;

1095:   options = options ? options : defaultoptions;
1096:   PetscOptionsValidKey(newname,&valid);
1097:   if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliased option %s",newname);
1098:   PetscOptionsValidKey(oldname,&valid);
1099:   if (!valid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid aliasee option %s",oldname);

1101:   n = options->Naliases;
1102:   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);

1104:   newname++; oldname++;
1105:   PetscStrlen(newname,&len);
1106:   options->aliases1[n] = (char*)malloc((len+1)*sizeof(char));
1107:   PetscStrcpy(options->aliases1[n],newname);
1108:   PetscStrlen(oldname,&len);
1109:   options->aliases2[n] = (char*)malloc((len+1)*sizeof(char));
1110:   PetscStrcpy(options->aliases2[n],oldname);
1111:   options->Naliases++;
1112:   return(0);
1113: }

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

1119:    Logically Collective

1121:    Input Parameters:
1122: +  options - options database, use NULL for the default global database
1123: .  name - name of option, this SHOULD have the - prepended
1124: -  value - the option value (not used for all options, so can be NULL)

1126:    Level: intermediate

1128:    Note:
1129:    This function can be called BEFORE PetscInitialize()

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

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

1138: .seealso: PetscOptionsInsert(), PetscOptionsClearValue()
1139: @*/
1140: PetscErrorCode PetscOptionsSetValue(PetscOptions options,const char name[],const char value[])
1141: {
1142:   return PetscOptionsSetValue_Private(options,name,value,NULL);
1143: }

1145: static PetscErrorCode PetscOptionsSetValue_Private(PetscOptions options,const char name[],const char value[],int *pos)
1146: {
1147:   size_t         len;
1148:   int            N,n,i;
1149:   char           **names;
1150:   char           fullname[MAXOPTNAME] = "";
1151:   PetscBool      flg;

1154:   if (!options) {
1155:     PetscOptionsCreateDefault();if (ierr) return ierr;
1156:     options = defaultoptions;
1157:   }

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

1161:   PetscOptionsSkipPrecedent(options,name,&flg);
1162:   if (flg) return 0;

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

1166:   if (options->prefixind > 0) {
1167:     strncpy(fullname,options->prefix,sizeof(fullname));
1168:     fullname[sizeof(fullname)-1] = 0;
1169:     strncat(fullname,name,sizeof(fullname)-strlen(fullname)-1);
1170:     fullname[sizeof(fullname)-1] = 0;
1171:     name = fullname;
1172:   }

1174:   /* check against aliases */
1175:   N = options->Naliases;
1176:   for (i=0; i<N; i++) {
1177:     int result = PetscOptNameCmp(options->aliases1[i],name);
1178:     if (!result) { name = options->aliases2[i]; break; }
1179:   }

1181:   /* slow search */
1182:   N = n = options->N;
1183:   names = options->names;
1184:   for (i=0; i<N; i++) {
1185:     int result = PetscOptNameCmp(names[i],name);
1186:     if (!result) {
1187:       n = i; goto setvalue;
1188:     } else if (result > 0) {
1189:       n = i; break;
1190:     }
1191:   }
1192:   if (N >= MAXOPTIONS) return PETSC_ERR_MEM;
1193:   /* shift remaining values up 1 */
1194:   for (i=N; i>n; i--) {
1195:     options->names[i]  = options->names[i-1];
1196:     options->values[i] = options->values[i-1];
1197:     options->used[i]   = options->used[i-1];
1198:   }
1199:   options->names[n]  = NULL;
1200:   options->values[n] = NULL;
1201:   options->used[n]   = PETSC_FALSE;
1202:   options->N++;

1204:   /* destroy hash table */
1205:   kh_destroy(HO,options->ht);
1206:   options->ht = NULL;

1208:   /* set new name */
1209:   len = strlen(name);
1210:   options->names[n] = (char*)malloc((len+1)*sizeof(char));
1211:   if (!options->names[n]) return PETSC_ERR_MEM;
1212:   strcpy(options->names[n],name);

1214: setvalue:
1215:   /* set new value */
1216:   if (options->values[n]) free(options->values[n]);
1217:   len = value ? strlen(value) : 0;
1218:   if (len) {
1219:     options->values[n] = (char*)malloc((len+1)*sizeof(char));
1220:     if (!options->values[n]) return PETSC_ERR_MEM;
1221:     strcpy(options->values[n],value);
1222:   } else {
1223:     options->values[n] = NULL;
1224:   }

1226:   /* handle -help so that it can be set from anywhere */
1227:   if (!PetscOptNameCmp(name,"help")) {
1228:     options->help = PETSC_TRUE;
1229:     PetscStrcasecmp(value, "intro", &options->help_intro);
1230:     options->used[n] = PETSC_TRUE;
1231:   }

1233:   if (PetscErrorHandlingInitialized) {
1234:     PetscOptionsMonitor(options,name,value);
1235:   }
1236:   if (pos) *pos = n;
1237:   return 0;
1238: }

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

1244:    Logically Collective

1246:    Input Parameter:
1247: +  options - options database, use NULL for the default global database
1248: -  name - name of option, this SHOULD have the - prepended

1250:    Level: intermediate

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

1257: .seealso: PetscOptionsInsert()
1258: @*/
1259: PetscErrorCode PetscOptionsClearValue(PetscOptions options,const char name[])
1260: {
1261:   int            N,n,i;
1262:   char           **names;

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

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

1272:   /* slow search */
1273:   N = n = options->N;
1274:   names = options->names;
1275:   for (i=0; i<N; i++) {
1276:     int result = PetscOptNameCmp(names[i],name);
1277:     if (!result) {
1278:       n = i; break;
1279:     } else if (result > 0) {
1280:       n = N; break;
1281:     }
1282:   }
1283:   if (n == N) return(0); /* it was not present */

1285:   /* remove name and value */
1286:   if (options->names[n])  free(options->names[n]);
1287:   if (options->values[n]) free(options->values[n]);
1288:   /* shift remaining values down 1 */
1289:   for (i=n; i<N-1; i++) {
1290:     options->names[i]  = options->names[i+1];
1291:     options->values[i] = options->values[i+1];
1292:     options->used[i]   = options->used[i+1];
1293:   }
1294:   options->N--;

1296:   /* destroy hash table */
1297:   kh_destroy(HO,options->ht);
1298:   options->ht = NULL;

1300:   PetscOptionsMonitor(options,name,NULL);
1301:   return(0);
1302: }

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

1307:    Not Collective

1309:    Input Parameters:
1310: +  options - options database, use NULL for the default global database
1311: .  pre - the string to prepend to the name or NULL, this SHOULD NOT have the "-" prepended
1312: -  name - name of option, this SHOULD have the "-" prepended

1314:    Output Parameters:
1315: +  value - the option value (optional, not used for all options)
1316: -  set - whether the option is set (optional)

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

1321:    Level: developer

1323: .seealso: PetscOptionsSetValue(), PetscOptionsClearValue()
1324: @*/
1325: PetscErrorCode PetscOptionsFindPair(PetscOptions options,const char pre[],const char name[],const char *value[],PetscBool *set)
1326: {
1327:   char           buf[MAXOPTNAME];
1328:   PetscBool      usehashtable = PETSC_TRUE;
1329:   PetscBool      matchnumbers = PETSC_TRUE;

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

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

1339:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1340:   if (pre && pre[0]) {
1341:     char *ptr = buf;
1342:     if (name[0] == '-') { *ptr++ = '-';  name++; }
1343:     PetscStrncpy(ptr,pre,buf+sizeof(buf)-ptr);
1344:     PetscStrlcat(buf,name,sizeof(buf));
1345:     name = buf;
1346:   }

1348:   if (PetscDefined(USE_DEBUG)) {
1349:     PetscBool valid;
1350:     char      key[MAXOPTNAME+1] = "-";
1351:     PetscStrncpy(key+1,name,sizeof(key)-1);
1352:     PetscOptionsValidKey(key,&valid);
1353:     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1354:   }

1356:   if (!options->ht && usehashtable) {
1357:     int i,ret;
1358:     khiter_t it;
1359:     khash_t(HO) *ht;
1360:     ht = kh_init(HO);
1361:     if (PetscUnlikely(!ht)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1362:     ret = kh_resize(HO,ht,options->N*2); /* twice the required size to reduce risk of collisions */
1363:     if (PetscUnlikely(ret)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1364:     for (i=0; i<options->N; i++) {
1365:       it = kh_put(HO,ht,options->names[i],&ret);
1366:       if (PetscUnlikely(ret != 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed");
1367:       kh_val(ht,it) = i;
1368:     }
1369:     options->ht = ht;
1370:   }

1372:   if (usehashtable)
1373:   { /* fast search */
1374:     khash_t(HO) *ht = options->ht;
1375:     khiter_t it = kh_get(HO,ht,name);
1376:     if (it != kh_end(ht)) {
1377:       int i = kh_val(ht,it);
1378:       options->used[i]  = PETSC_TRUE;
1379:       if (value) *value = options->values[i];
1380:       if (set)   *set   = PETSC_TRUE;
1381:       return(0);
1382:     }
1383:   } else
1384:   { /* slow search */
1385:     int i, N = options->N;
1386:     for (i=0; i<N; i++) {
1387:       int result = PetscOptNameCmp(options->names[i],name);
1388:       if (!result) {
1389:         options->used[i]  = PETSC_TRUE;
1390:         if (value) *value = options->values[i];
1391:         if (set)   *set   = PETSC_TRUE;
1392:         return(0);
1393:       } else if (result > 0) {
1394:         break;
1395:       }
1396:     }
1397:   }

1399:   /*
1400:    The following block slows down all lookups in the most frequent path (most lookups are unsuccessful).
1401:    Maybe this special lookup mode should be enabled on request with a push/pop API.
1402:    The feature of matching _%d_ used sparingly in the codebase.
1403:    */
1404:   if (matchnumbers) {
1405:     int i,j,cnt = 0,locs[16],loce[16];
1406:     /* determine the location and number of all _%d_ in the key */
1407:     for (i=0; name[i]; i++) {
1408:       if (name[i] == '_') {
1409:         for (j=i+1; name[j]; j++) {
1410:           if (name[j] >= '0' && name[j] <= '9') continue;
1411:           if (name[j] == '_' && j > i+1) { /* found a number */
1412:             locs[cnt]   = i+1;
1413:             loce[cnt++] = j+1;
1414:           }
1415:           i = j-1;
1416:           break;
1417:         }
1418:       }
1419:     }
1420:     for (i=0; i<cnt; i++) {
1421:       PetscBool found;
1422:       char      opt[MAXOPTNAME+1] = "-", tmp[MAXOPTNAME];
1423:       PetscStrncpy(tmp,name,PetscMin((size_t)(locs[i]+1),sizeof(tmp)));
1424:       PetscStrlcat(opt,tmp,sizeof(opt));
1425:       PetscStrlcat(opt,name+loce[i],sizeof(opt));
1426:       PetscOptionsFindPair(options,NULL,opt,value,&found);
1427:       if (found) {if (set) *set = PETSC_TRUE; return(0);}
1428:     }
1429:   }

1431:   if (set) *set = PETSC_FALSE;
1432:   return(0);
1433: }

1435: /* Check whether any option begins with pre+name */
1436: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[],const char *value[],PetscBool *set)
1437: {
1438:   char           buf[MAXOPTNAME];
1439:   int            numCnt = 0, locs[16],loce[16];

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

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

1449:   /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */
1450:   if (pre && pre[0]) {
1451:     char *ptr = buf;
1452:     if (name[0] == '-') { *ptr++ = '-';  name++; }
1453:     PetscStrncpy(ptr,pre,sizeof(buf)+(size_t)(ptr-buf));
1454:     PetscStrlcat(buf,name,sizeof(buf));
1455:     name = buf;
1456:   }

1458:   if (PetscDefined(USE_DEBUG)) {
1459:     PetscBool valid;
1460:     char      key[MAXOPTNAME+1] = "-";
1461:     PetscStrncpy(key+1,name,sizeof(key)-1);
1462:     PetscOptionsValidKey(key,&valid);
1463:     if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name);
1464:   }

1466:   /* determine the location and number of all _%d_ in the key */
1467:   {
1468:     int i,j;
1469:     for (i=0; name[i]; i++) {
1470:       if (name[i] == '_') {
1471:         for (j=i+1; name[j]; j++) {
1472:           if (name[j] >= '0' && name[j] <= '9') continue;
1473:           if (name[j] == '_' && j > i+1) { /* found a number */
1474:             locs[numCnt]   = i+1;
1475:             loce[numCnt++] = j+1;
1476:           }
1477:           i = j-1;
1478:           break;
1479:         }
1480:       }
1481:     }
1482:   }

1484:   { /* slow search */
1485:     int       c, i;
1486:     size_t    len;
1487:     PetscBool match;

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

1492:       if (c < 0) {
1493:         PetscStrcpy(opt,name);
1494:       } else {
1495:         PetscStrncpy(tmp,name,PetscMin((size_t)(locs[c]+1),sizeof(tmp)));
1496:         PetscStrlcat(opt,tmp,sizeof(opt));
1497:         PetscStrlcat(opt,name+loce[c],sizeof(opt));
1498:       }
1499:       PetscStrlen(opt,&len);
1500:       for (i=0; i<options->N; i++) {
1501:         PetscStrncmp(options->names[i],opt,len,&match);
1502:         if (match) {
1503:           options->used[i]  = PETSC_TRUE;
1504:           if (value) *value = options->values[i];
1505:           if (set)   *set   = PETSC_TRUE;
1506:           return(0);
1507:         }
1508:       }
1509:     }
1510:   }

1512:   if (set) *set = PETSC_FALSE;
1513:   return(0);
1514: }

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

1519:    Not Collective

1521:    Input Parameters:
1522: +  options - options database, use NULL for default global database
1523: .  pre - the option prefix (may be NULL)
1524: .  name - the option name one is seeking
1525: -  mess - error message (may be NULL)

1527:    Level: advanced

1529: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(),
1530:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1531:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1532:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1533:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1534:           PetscOptionsFList(), PetscOptionsEList()
1535: @*/
1536: PetscErrorCode PetscOptionsReject(PetscOptions options,const char pre[],const char name[],const char mess[])
1537: {
1539:   PetscBool      flag = PETSC_FALSE;

1542:   PetscOptionsHasName(options,pre,name,&flag);
1543:   if (flag) {
1544:     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);
1545:     else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s",pre?pre:"",name+1);
1546:   }
1547:   return(0);
1548: }

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

1553:    Not Collective

1555:    Input Parameters:
1556: .  options - options database, use NULL for default global database

1558:    Output Parameters:
1559: .  set - PETSC_TRUE if found else PETSC_FALSE.

1561:    Level: advanced

1563: .seealso: PetscOptionsHasName()
1564: @*/
1565: PetscErrorCode PetscOptionsHasHelp(PetscOptions options,PetscBool *set)
1566: {
1569:   options = options ? options : defaultoptions;
1570:   *set = options->help;
1571:   return(0);
1572: }

1574: PetscErrorCode PetscOptionsHasHelpIntro_Internal(PetscOptions options,PetscBool *set)
1575: {
1578:   options = options ? options : defaultoptions;
1579:   *set = options->help_intro;
1580:   return(0);
1581: }

1583: /*@C
1584:    PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or boolean, even
1585:                       its value is set to false.

1587:    Not Collective

1589:    Input Parameters:
1590: +  options - options database, use NULL for default global database
1591: .  pre - string to prepend to the name or NULL
1592: -  name - the option one is seeking

1594:    Output Parameters:
1595: .  set - PETSC_TRUE if found else PETSC_FALSE.

1597:    Level: beginner

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

1602: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1603:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1604:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1605:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1606:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1607:           PetscOptionsFList(), PetscOptionsEList()
1608: @*/
1609: PetscErrorCode PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool *set)
1610: {
1611:   const char     *value;
1613:   PetscBool      flag;

1616:   PetscOptionsFindPair(options,pre,name,&value,&flag);
1617:   if (set) *set = flag;
1618:   return(0);
1619: }

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

1624:    Not Collective

1626:    Input Parameter:
1627: .  options - the options database, use NULL for the default global database

1629:    Output Parameter:
1630: .  copts - pointer where string pointer is stored

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

1636:    Level: advanced

1638: .seealso: PetscOptionsAllUsed(), PetscOptionsView(), PetscOptionsPush(), PetscOptionsPop()
1639: @*/
1640: PetscErrorCode PetscOptionsGetAll(PetscOptions options,char *copts[])
1641: {
1643:   PetscInt       i;
1644:   size_t         len = 1,lent = 0;
1645:   char           *coptions = NULL;

1649:   options = options ? options : defaultoptions;
1650:   /* count the length of the required string */
1651:   for (i=0; i<options->N; i++) {
1652:     PetscStrlen(options->names[i],&lent);
1653:     len += 2 + lent;
1654:     if (options->values[i]) {
1655:       PetscStrlen(options->values[i],&lent);
1656:       len += 1 + lent;
1657:     }
1658:   }
1659:   PetscMalloc1(len,&coptions);
1660:   coptions[0] = 0;
1661:   for (i=0; i<options->N; i++) {
1662:     PetscStrcat(coptions,"-");
1663:     PetscStrcat(coptions,options->names[i]);
1664:     PetscStrcat(coptions," ");
1665:     if (options->values[i]) {
1666:       PetscStrcat(coptions,options->values[i]);
1667:       PetscStrcat(coptions," ");
1668:     }
1669:   }
1670:   *copts = coptions;
1671:   return(0);
1672: }

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

1677:    Not Collective

1679:    Input Parameter:
1680: +  options - options database, use NULL for default global database
1681: -  name - string name of option

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

1686:    Level: advanced

1688:    Notes:
1689:    The value returned may be different on each process and depends on which options have been processed
1690:    on the given process

1692: .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed()
1693: @*/
1694: PetscErrorCode PetscOptionsUsed(PetscOptions options,const char *name,PetscBool *used)
1695: {
1696:   PetscInt       i;

1702:   options = options ? options : defaultoptions;
1703:   *used = PETSC_FALSE;
1704:   for (i=0; i<options->N; i++) {
1705:     PetscStrcasecmp(options->names[i],name,used);
1706:     if (*used) {
1707:       *used = options->used[i];
1708:       break;
1709:     }
1710:   }
1711:   return(0);
1712: }

1714: /*@
1715:    PetscOptionsAllUsed - Returns a count of the number of options in the
1716:    database that have never been selected.

1718:    Not Collective

1720:    Input Parameter:
1721: .  options - options database, use NULL for default global database

1723:    Output Parameter:
1724: .  N - count of options not used

1726:    Level: advanced

1728:    Notes:
1729:    The value returned may be different on each process and depends on which options have been processed
1730:    on the given process

1732: .seealso: PetscOptionsView()
1733: @*/
1734: PetscErrorCode PetscOptionsAllUsed(PetscOptions options,PetscInt *N)
1735: {
1736:   PetscInt     i,n = 0;

1740:   options = options ? options : defaultoptions;
1741:   for (i=0; i<options->N; i++) {
1742:     if (!options->used[i]) n++;
1743:   }
1744:   *N = n;
1745:   return(0);
1746: }

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

1751:    Not Collective

1753:    Input Parameter:
1754: .  options - options database; use NULL for default global database

1756:    Options Database Key:
1757: .  -options_left - activates PetscOptionsAllUsed() within PetscFinalize()

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

1765:    Level: advanced

1767: .seealso: PetscOptionsAllUsed()
1768: @*/
1769: PetscErrorCode PetscOptionsLeft(PetscOptions options)
1770: {
1772:   PetscInt       i;
1773:   PetscInt       cnt = 0;
1774:   PetscOptions   toptions;

1777:   toptions = options ? options : defaultoptions;
1778:   for (i=0; i<toptions->N; i++) {
1779:     if (!toptions->used[i]) {
1780:       if (toptions->values[i]) {
1781:         PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",toptions->names[i],toptions->values[i]);
1782:       } else {
1783:         PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",toptions->names[i]);
1784:       }
1785:     }
1786:   }
1787:   if (!options) {
1788:     toptions = defaultoptions;
1789:     while (toptions->previous) {
1790:       cnt++;
1791:       toptions = toptions->previous;
1792:     }
1793:     if (cnt) {
1794:       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);
1795:     }
1796:   }
1797:   return(0);
1798: }

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

1803:    Not Collective

1805:    Input Parameter:
1806: .  options - options database, use NULL for default global database

1808:    Output Parameter:
1809: +  N - count of options not used
1810: .  names - names of options not used
1811: -  values - values of options not used

1813:    Level: advanced

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

1820: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft()
1821: @*/
1822: PetscErrorCode PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[])
1823: {
1825:   PetscInt       i,n;

1831:   options = options ? options : defaultoptions;

1833:   /* The number of unused PETSc options */
1834:   n = 0;
1835:   for (i=0; i<options->N; i++) {
1836:     if (!options->used[i]) n++;
1837:   }
1838:   if (N) { *N = n; }
1839:   if (names)  { PetscMalloc1(n,names); }
1840:   if (values) { PetscMalloc1(n,values); }

1842:   n = 0;
1843:   if (names || values) {
1844:     for (i=0; i<options->N; i++) {
1845:       if (!options->used[i]) {
1846:         if (names)  (*names)[n]  = options->names[i];
1847:         if (values) (*values)[n] = options->values[i];
1848:         n++;
1849:       }
1850:     }
1851:   }
1852:   return(0);
1853: }

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

1858:    Not Collective

1860:    Input Parameter:
1861: +  options - options database, use NULL for default global database
1862: .  names - names of options not used
1863: -  values - values of options not used

1865:    Level: advanced

1867: .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet()
1868: @*/
1869: PetscErrorCode PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[])
1870: {

1877:   if (N) { *N = 0; }
1878:   if (names)  { PetscFree(*names); }
1879:   if (values) { PetscFree(*values); }
1880:   return(0);
1881: }

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

1886:    Logically Collective on ctx

1888:    Input Parameters:
1889: +  name  - option name string
1890: .  value - option value string
1891: -  ctx - an ASCII viewer or NULL

1893:    Level: intermediate

1895:    Notes:
1896:      If ctx=NULL, PetscPrintf() is used.
1897:      The first MPI rank in the PetscViewer viewer actually prints the values, other
1898:      processes may have different values set

1900: .seealso: PetscOptionsMonitorSet()
1901: @*/
1902: PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx)
1903: {

1907:   if (ctx) {
1908:     PetscViewer viewer = (PetscViewer)ctx;
1909:     if (!value) {
1910:       PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);
1911:     } else if (!value[0]) {
1912:       PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);
1913:     } else {
1914:       PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);
1915:     }
1916:   } else {
1917:     MPI_Comm comm = PETSC_COMM_WORLD;
1918:     if (!value) {
1919:       PetscPrintf(comm,"Removing option: %s\n",name,value);
1920:     } else if (!value[0]) {
1921:       PetscPrintf(comm,"Setting option: %s (no value)\n",name);
1922:     } else {
1923:       PetscPrintf(comm,"Setting option: %s = %s\n",name,value);
1924:     }
1925:   }
1926:   return(0);
1927: }

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

1933:    Not Collective

1935:    Input Parameters:
1936: +  monitor - pointer to function (if this is NULL, it turns off monitoring
1937: .  mctx    - [optional] context for private data for the
1938:              monitor routine (use NULL if no context is desired)
1939: -  monitordestroy - [optional] routine that frees monitor context
1940:           (may be NULL)

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

1945: +  name - option name string
1946: .  value - option value string
1947: -  mctx  - optional monitoring context, as set by PetscOptionsMonitorSet()

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

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

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

1961:    Level: intermediate

1963: .seealso: PetscOptionsMonitorDefault(), PetscInitialize()
1964: @*/
1965: PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1966: {
1967:   PetscOptions options = defaultoptions;

1970:   if (options->monitorCancel) return(0);
1971:   if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set");
1972:   options->monitor[options->numbermonitors]          = monitor;
1973:   options->monitordestroy[options->numbermonitors]   = monitordestroy;
1974:   options->monitorcontext[options->numbermonitors++] = (void*)mctx;
1975:   return(0);
1976: }

1978: /*
1979:    PetscOptionsStringToBool - Converts string to PetscBool, handles cases like "yes", "no", "true", "false", "0", "1", "off", "on".
1980:      Empty string is considered as true.
1981: */
1982: PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a)
1983: {
1984:   PetscBool      istrue,isfalse;
1985:   size_t         len;

1989:   /* PetscStrlen() returns 0 for NULL or "" */
1990:   PetscStrlen(value,&len);
1991:   if (!len)  {*a = PETSC_TRUE; return(0);}
1992:   PetscStrcasecmp(value,"TRUE",&istrue);
1993:   if (istrue) {*a = PETSC_TRUE; return(0);}
1994:   PetscStrcasecmp(value,"YES",&istrue);
1995:   if (istrue) {*a = PETSC_TRUE; return(0);}
1996:   PetscStrcasecmp(value,"1",&istrue);
1997:   if (istrue) {*a = PETSC_TRUE; return(0);}
1998:   PetscStrcasecmp(value,"on",&istrue);
1999:   if (istrue) {*a = PETSC_TRUE; return(0);}
2000:   PetscStrcasecmp(value,"FALSE",&isfalse);
2001:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2002:   PetscStrcasecmp(value,"NO",&isfalse);
2003:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2004:   PetscStrcasecmp(value,"0",&isfalse);
2005:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2006:   PetscStrcasecmp(value,"off",&isfalse);
2007:   if (isfalse) {*a = PETSC_FALSE; return(0);}
2008:   SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value);
2009: }

2011: /*
2012:    PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide"
2013: */
2014: PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a)
2015: {
2017:   size_t         len;
2018:   PetscBool      decide,tdefault,mouse;

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

2024:   PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);
2025:   if (!tdefault) {
2026:     PetscStrcasecmp(name,"DEFAULT",&tdefault);
2027:   }
2028:   PetscStrcasecmp(name,"PETSC_DECIDE",&decide);
2029:   if (!decide) {
2030:     PetscStrcasecmp(name,"DECIDE",&decide);
2031:   }
2032:   PetscStrcasecmp(name,"mouse",&mouse);

2034:   if (tdefault)    *a = PETSC_DEFAULT;
2035:   else if (decide) *a = PETSC_DECIDE;
2036:   else if (mouse)  *a = -1;
2037:   else {
2038:     char *endptr;
2039:     long strtolval;

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

2044: #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL)
2045:     (void) strtolval;
2046:     *a = atoll(name);
2047: #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64)
2048:     (void) strtolval;
2049:     *a = _atoi64(name);
2050: #else
2051:     *a = (PetscInt)strtolval;
2052: #endif
2053:   }
2054:   return(0);
2055: }

2057: #if defined(PETSC_USE_REAL___FLOAT128)
2058: #include <quadmath.h>
2059: #endif

2061: static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr)
2062: {
2064: #if defined(PETSC_USE_REAL___FLOAT128)
2065:   *a = strtoflt128(name,endptr);
2066: #else
2067:   *a = (PetscReal)strtod(name,endptr);
2068: #endif
2069:   return(0);
2070: }

2072: static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary)
2073: {
2074:   PetscBool      hasi = PETSC_FALSE;
2075:   char           *ptr;
2076:   PetscReal      strtoval;

2080:   PetscStrtod(name,&strtoval,&ptr);
2081:   if (ptr == name) {
2082:     strtoval = 1.;
2083:     hasi = PETSC_TRUE;
2084:     if (name[0] == 'i') {
2085:       ptr++;
2086:     } else if (name[0] == '+' && name[1] == 'i') {
2087:       ptr += 2;
2088:     } else if (name[0] == '-' && name[1] == 'i') {
2089:       strtoval = -1.;
2090:       ptr += 2;
2091:     }
2092:   } else if (*ptr == 'i') {
2093:     hasi = PETSC_TRUE;
2094:     ptr++;
2095:   }
2096:   *endptr = ptr;
2097:   *isImaginary = hasi;
2098:   if (hasi) {
2099: #if !defined(PETSC_USE_COMPLEX)
2100:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name);
2101: #else
2102:     *a = PetscCMPLX(0.,strtoval);
2103: #endif
2104:   } else {
2105:     *a = strtoval;
2106:   }
2107:   return(0);
2108: }

2110: /*
2111:    Converts a string to PetscReal value. Handles special cases like "default" and "decide"
2112: */
2113: PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a)
2114: {
2115:   size_t         len;
2116:   PetscBool      match;
2117:   char           *endptr;

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

2124:   PetscStrcasecmp(name,"PETSC_DEFAULT",&match);
2125:   if (!match) {
2126:     PetscStrcasecmp(name,"DEFAULT",&match);
2127:   }
2128:   if (match) {*a = PETSC_DEFAULT; return(0);}

2130:   PetscStrcasecmp(name,"PETSC_DECIDE",&match);
2131:   if (!match) {
2132:     PetscStrcasecmp(name,"DECIDE",&match);
2133:   }
2134:   if (match) {*a = PETSC_DECIDE; return(0);}

2136:   PetscStrtod(name,a,&endptr);
2137:   if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name);
2138:   return(0);
2139: }

2141: PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a)
2142: {
2143:   PetscBool      imag1;
2144:   size_t         len;
2145:   PetscScalar    val = 0.;
2146:   char           *ptr = NULL;

2150:   PetscStrlen(name,&len);
2151:   if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value");
2152:   PetscStrtoz(name,&val,&ptr,&imag1);
2153: #if defined(PETSC_USE_COMPLEX)
2154:   if ((size_t) (ptr - name) < len) {
2155:     PetscBool   imag2;
2156:     PetscScalar val2;

2158:     PetscStrtoz(ptr,&val2,&ptr,&imag2);
2159:     if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name);
2160:     val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2));
2161:   }
2162: #endif
2163:   if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name);
2164:   *a = val;
2165:   return(0);
2166: }

2168: /*@C
2169:    PetscOptionsGetBool - Gets the Logical (true or false) value for a particular
2170:             option in the database.

2172:    Not Collective

2174:    Input Parameters:
2175: +  options - options database, use NULL for default global database
2176: .  pre - the string to prepend to the name or NULL
2177: -  name - the option one is seeking

2179:    Output Parameter:
2180: +  ivalue - the logical value to return
2181: -  set - PETSC_TRUE  if found, else PETSC_FALSE

2183:    Level: beginner

2185:    Notes:
2186:        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2187:        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE

2189:       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
2190:      is equivalent to -requested_bool true

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

2195: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
2196:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(),
2197:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2198:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2199:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2200:           PetscOptionsFList(), PetscOptionsEList()
2201: @*/
2202: PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set)
2203: {
2204:   const char     *value;
2205:   PetscBool      flag;

2211:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2212:   if (flag) {
2213:     if (set) *set = PETSC_TRUE;
2214:     PetscOptionsStringToBool(value, &flag);
2215:     if (ivalue) *ivalue = flag;
2216:   } else {
2217:     if (set) *set = PETSC_FALSE;
2218:   }
2219:   return(0);
2220: }

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

2225:    Not Collective

2227:    Input Parameters:
2228: +  options - options database, use NULL for default global database
2229: .  pre - the string to prepend to the name or NULL
2230: .  opt - option name
2231: .  list - the possible choices (one of these must be selected, anything else is invalid)
2232: -  ntext - number of choices

2234:    Output Parameter:
2235: +  value - the index of the value to return (defaults to zero if the option name is given but no choice is listed)
2236: -  set - PETSC_TRUE if found, else PETSC_FALSE

2238:    Level: intermediate

2240:    Notes:
2241:     If the user does not supply the option value 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:    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()

2246: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2247:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2248:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2249:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2250:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2251:           PetscOptionsFList(), PetscOptionsEList()
2252: @*/
2253: PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set)
2254: {
2256:   size_t         alen,len = 0, tlen = 0;
2257:   char           *svalue;
2258:   PetscBool      aset,flg = PETSC_FALSE;
2259:   PetscInt       i;

2263:   for (i=0; i<ntext; i++) {
2264:     PetscStrlen(list[i],&alen);
2265:     if (alen > len) len = alen;
2266:     tlen += len + 1;
2267:   }
2268:   len += 5; /* a little extra space for user mistypes */
2269:   PetscMalloc1(len,&svalue);
2270:   PetscOptionsGetString(options,pre,opt,svalue,len,&aset);
2271:   if (aset) {
2272:     PetscEListFind(ntext,list,svalue,value,&flg);
2273:     if (!flg) {
2274:       char *avail,*pavl;

2276:       PetscMalloc1(tlen,&avail);
2277:       pavl = avail;
2278:       for (i=0; i<ntext; i++) {
2279:         PetscStrlen(list[i],&alen);
2280:         PetscStrcpy(pavl,list[i]);
2281:         pavl += alen;
2282:         PetscStrcpy(pavl," ");
2283:         pavl += 1;
2284:       }
2285:       PetscStrtolower(avail);
2286:       SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s. Available options: %s",svalue,pre ? pre : "",opt+1,avail);
2287:     }
2288:     if (set) *set = PETSC_TRUE;
2289:   } else if (set) *set = PETSC_FALSE;
2290:   PetscFree(svalue);
2291:   return(0);
2292: }

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

2297:    Not Collective

2299:    Input Parameters:
2300: +  options - options database, use NULL for default global database
2301: .  pre - option prefix or NULL
2302: .  opt - option name
2303: .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2304: -  defaultv - the default (current) value

2306:    Output Parameter:
2307: +  value - the  value to return
2308: -  set - PETSC_TRUE if found, else PETSC_FALSE

2310:    Level: beginner

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

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

2318: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2319:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2320:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
2321:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2322:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2323:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2324:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2325: @*/
2326: PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set)
2327: {
2329:   PetscInt       ntext = 0,tval;
2330:   PetscBool      fset;

2334:   while (list[ntext++]) {
2335:     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
2336:   }
2337:   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
2338:   ntext -= 3;
2339:   PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);
2340:   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
2341:   if (fset) *value = (PetscEnum)tval;
2342:   if (set) *set = fset;
2343:   return(0);
2344: }

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

2349:    Not Collective

2351:    Input Parameters:
2352: +  options - options database, use NULL for default global database
2353: .  pre - the string to prepend to the name or NULL
2354: -  name - the option one is seeking

2356:    Output Parameter:
2357: +  ivalue - the integer value to return
2358: -  set - PETSC_TRUE if found, else PETSC_FALSE

2360:    Level: beginner

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

2366: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
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()
2373: @*/
2374: PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set)
2375: {
2376:   const char     *value;
2378:   PetscBool      flag;

2383:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2384:   if (flag) {
2385:     if (!value) {
2386:       if (set) *set = PETSC_FALSE;
2387:     } else {
2388:       if (set) *set = PETSC_TRUE;
2389:       PetscOptionsStringToInt(value,ivalue);
2390:     }
2391:   } else {
2392:     if (set) *set = PETSC_FALSE;
2393:   }
2394:   return(0);
2395: }

2397: /*@C
2398:    PetscOptionsGetReal - Gets the double precision value for a particular
2399:    option in the database.

2401:    Not Collective

2403:    Input Parameters:
2404: +  options - options database, use NULL for default global database
2405: .  pre - string to prepend to each name or NULL
2406: -  name - the option one is seeking

2408:    Output Parameter:
2409: +  dvalue - the double value to return
2410: -  set - PETSC_TRUE if found, PETSC_FALSE if not found

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

2416:    Level: beginner

2418: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2419:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(),
2420:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2421:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2422:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2423:           PetscOptionsFList(), PetscOptionsEList()
2424: @*/
2425: PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set)
2426: {
2427:   const char     *value;
2428:   PetscBool      flag;

2434:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2435:   if (flag) {
2436:     if (!value) {
2437:       if (set) *set = PETSC_FALSE;
2438:     } else {
2439:       if (set) *set = PETSC_TRUE;
2440:       PetscOptionsStringToReal(value,dvalue);
2441:     }
2442:   } else {
2443:     if (set) *set = PETSC_FALSE;
2444:   }
2445:   return(0);
2446: }

2448: /*@C
2449:    PetscOptionsGetScalar - Gets the scalar value for a particular
2450:    option in the database.

2452:    Not Collective

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

2459:    Output Parameter:
2460: +  dvalue - the double value to return
2461: -  set - PETSC_TRUE if found, else PETSC_FALSE

2463:    Level: beginner

2465:    Usage:
2466:    A complex number 2+3i must be specified with NO spaces

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

2472: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2473:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2474:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2475:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2476:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2477:           PetscOptionsFList(), PetscOptionsEList()
2478: @*/
2479: PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set)
2480: {
2481:   const char     *value;
2482:   PetscBool      flag;

2488:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2489:   if (flag) {
2490:     if (!value) {
2491:       if (set) *set = PETSC_FALSE;
2492:     } else {
2493: #if !defined(PETSC_USE_COMPLEX)
2494:       PetscOptionsStringToReal(value,dvalue);
2495: #else
2496:       PetscOptionsStringToScalar(value,dvalue);
2497: #endif
2498:       if (set) *set = PETSC_TRUE;
2499:     }
2500:   } else { /* flag */
2501:     if (set) *set = PETSC_FALSE;
2502:   }
2503:   return(0);
2504: }

2506: /*@C
2507:    PetscOptionsGetString - Gets the string value for a particular option in
2508:    the database.

2510:    Not Collective

2512:    Input Parameters:
2513: +  options - options database, use NULL for default global database
2514: .  pre - string to prepend to name or NULL
2515: .  name - the option one is seeking
2516: -  len - maximum length of the string including null termination

2518:    Output Parameters:
2519: +  string - location to copy string
2520: -  set - PETSC_TRUE if found, else PETSC_FALSE

2522:    Level: beginner

2524:    Fortran Note:
2525:    The Fortran interface is slightly different from the C/C++
2526:    interface (len is not used).  Sample usage in Fortran follows
2527: .vb
2528:       character *20    string
2529:       PetscErrorCode   ierr
2530:       PetscBool        set
2531:       call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr)
2532: .ve

2534:    Notes:
2535:     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

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

2540:     Note:
2541:       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).

2543: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2544:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2545:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2546:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2547:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2548:           PetscOptionsFList(), PetscOptionsEList()
2549: @*/
2550: PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set)
2551: {
2552:   const char     *value;
2553:   PetscBool      flag;

2559:   PetscOptionsFindPair(options,pre,name,&value,&flag);
2560:   if (!flag) {
2561:     if (set) *set = PETSC_FALSE;
2562:   } else {
2563:     if (set) *set = PETSC_TRUE;
2564:     if (value) {
2565:       PetscStrncpy(string,value,len);
2566:     } else {
2567:       PetscArrayzero(string,len);
2568:     }
2569:   }
2570:   return(0);
2571: }

2573: char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[])
2574: {
2575:   const char     *value;
2576:   PetscBool      flag;

2580:   PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(NULL);
2581:   if (flag) PetscFunctionReturn((char*)value);
2582:   else PetscFunctionReturn(NULL);
2583: }

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

2590:    Not Collective

2592:    Input Parameters:
2593: +  options - options database, use NULL for default global database
2594: .  pre - string to prepend to each name or NULL
2595: .  name - the option one is seeking
2596: -  nmax - maximum number of values to retrieve

2598:    Output Parameter:
2599: +  dvalue - the integer values to return
2600: .  nmax - actual number of values retreived
2601: -  set - PETSC_TRUE if found, else PETSC_FALSE

2603:    Level: beginner

2605:    Notes:
2606:        TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE
2607:        FALSE, false, NO, no, and 0 all translate to PETSC_FALSE

2609: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2610:           PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2611:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2612:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2613:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2614:           PetscOptionsFList(), PetscOptionsEList()
2615: @*/
2616: PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set)
2617: {
2618:   const char     *svalue;
2619:   char           *value;
2621:   PetscInt       n = 0;
2622:   PetscBool      flag;
2623:   PetscToken     token;


2630:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2631:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2632:   if (set) *set = PETSC_TRUE;
2633:   PetscTokenCreate(svalue,',',&token);
2634:   PetscTokenFind(token,&value);
2635:   while (value && n < *nmax) {
2636:     PetscOptionsStringToBool(value,dvalue);
2637:     PetscTokenFind(token,&value);
2638:     dvalue++;
2639:     n++;
2640:   }
2641:   PetscTokenDestroy(&token);
2642:   *nmax = n;
2643:   return(0);
2644: }

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

2649:    Not Collective

2651:    Input Parameters:
2652: +  options - options database, use NULL for default global database
2653: .  pre - option prefix or NULL
2654: .  name - option name
2655: .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
2656: -  nmax - maximum number of values to retrieve

2658:    Output Parameters:
2659: +  ivalue - the  enum values to return
2660: .  nmax - actual number of values retreived
2661: -  set - PETSC_TRUE if found, else PETSC_FALSE

2663:    Level: beginner

2665:    Notes:
2666:    The array must be passed as a comma separated list.

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

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

2672: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
2673:           PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
2674:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(),
2675:           PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(),
2676:           PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2677:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum()
2678: @*/
2679: PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set)
2680: {
2681:   const char     *svalue;
2682:   char           *value;
2683:   PetscInt       n = 0;
2684:   PetscEnum      evalue;
2685:   PetscBool      flag;
2686:   PetscToken     token;


2695:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2696:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2697:   if (set) *set = PETSC_TRUE;
2698:   PetscTokenCreate(svalue,',',&token);
2699:   PetscTokenFind(token,&value);
2700:   while (value && n < *nmax) {
2701:     PetscEnumFind(list,value,&evalue,&flag);
2702:     if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1);
2703:     ivalue[n++] = evalue;
2704:     PetscTokenFind(token,&value);
2705:   }
2706:   PetscTokenDestroy(&token);
2707:   *nmax = n;
2708:   return(0);
2709: }

2711: /*@C
2712:    PetscOptionsGetIntArray - Gets an array of integer values for a particular
2713:    option in the database.

2715:    Not Collective

2717:    Input Parameters:
2718: +  options - options database, use NULL for default global database
2719: .  pre - string to prepend to each name or NULL
2720: .  name - the option one is seeking
2721: -  nmax - maximum number of values to retrieve

2723:    Output Parameter:
2724: +  ivalue - the integer values to return
2725: .  nmax - actual number of values retreived
2726: -  set - PETSC_TRUE if found, else PETSC_FALSE

2728:    Level: beginner

2730:    Notes:
2731:    The array can be passed as
2732:    a comma separated list:                                 0,1,2,3,4,5,6,7
2733:    a range (start-end+1):                                  0-8
2734:    a range with given increment (start-end+1:inc):         0-7:2
2735:    a combination of values and ranges separated by commas: 0,1-8,8-15:2

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

2739: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2740:           PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2741:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2742:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2743:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2744:           PetscOptionsFList(), PetscOptionsEList()
2745: @*/
2746: PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set)
2747: {
2748:   const char     *svalue;
2749:   char           *value;
2751:   PetscInt       n = 0,i,j,start,end,inc,nvalues;
2752:   size_t         len;
2753:   PetscBool      flag,foundrange;
2754:   PetscToken     token;


2761:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2762:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2763:   if (set) *set = PETSC_TRUE;
2764:   PetscTokenCreate(svalue,',',&token);
2765:   PetscTokenFind(token,&value);
2766:   while (value && n < *nmax) {
2767:     /* look for form  d-D where d and D are integers */
2768:     foundrange = PETSC_FALSE;
2769:     PetscStrlen(value,&len);
2770:     if (value[0] == '-') i=2;
2771:     else i=1;
2772:     for (;i<(int)len; i++) {
2773:       if (value[i] == '-') {
2774:         if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
2775:         value[i] = 0;

2777:         PetscOptionsStringToInt(value,&start);
2778:         inc  = 1;
2779:         j    = i+1;
2780:         for (;j<(int)len; j++) {
2781:           if (value[j] == ':') {
2782:             value[j] = 0;

2784:             PetscOptionsStringToInt(value+j+1,&inc);
2785:             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);
2786:             break;
2787:           }
2788:         }
2789:         PetscOptionsStringToInt(value+i+1,&end);
2790:         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);
2791:         nvalues = (end-start)/inc + (end-start)%inc;
2792:         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);
2793:         for (;start<end; start+=inc) {
2794:           *ivalue = start; ivalue++;n++;
2795:         }
2796:         foundrange = PETSC_TRUE;
2797:         break;
2798:       }
2799:     }
2800:     if (!foundrange) {
2801:       PetscOptionsStringToInt(value,ivalue);
2802:       ivalue++;
2803:       n++;
2804:     }
2805:     PetscTokenFind(token,&value);
2806:   }
2807:   PetscTokenDestroy(&token);
2808:   *nmax = n;
2809:   return(0);
2810: }

2812: /*@C
2813:    PetscOptionsGetRealArray - Gets an array of double precision values for a
2814:    particular option in the database.  The values must be separated with
2815:    commas with no intervening spaces.

2817:    Not Collective

2819:    Input Parameters:
2820: +  options - options database, use NULL for default global database
2821: .  pre - string to prepend to each name or NULL
2822: .  name - the option one is seeking
2823: -  nmax - maximum number of values to retrieve

2825:    Output Parameters:
2826: +  dvalue - the double values to return
2827: .  nmax - actual number of values retreived
2828: -  set - PETSC_TRUE if found, else PETSC_FALSE

2830:    Level: beginner

2832: .seealso: PetscOptionsGetInt(), PetscOptionsHasName(),
2833:           PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(),
2834:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2835:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2836:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2837:           PetscOptionsFList(), PetscOptionsEList()
2838: @*/
2839: PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set)
2840: {
2841:   const char     *svalue;
2842:   char           *value;
2844:   PetscInt       n = 0;
2845:   PetscBool      flag;
2846:   PetscToken     token;


2853:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2854:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2855:   if (set) *set = PETSC_TRUE;
2856:   PetscTokenCreate(svalue,',',&token);
2857:   PetscTokenFind(token,&value);
2858:   while (value && n < *nmax) {
2859:     PetscOptionsStringToReal(value,dvalue++);
2860:     PetscTokenFind(token,&value);
2861:     n++;
2862:   }
2863:   PetscTokenDestroy(&token);
2864:   *nmax = n;
2865:   return(0);
2866: }

2868: /*@C
2869:    PetscOptionsGetScalarArray - Gets an array of scalars for a
2870:    particular option in the database.  The values must be separated with
2871:    commas with no intervening spaces.

2873:    Not Collective

2875:    Input Parameters:
2876: +  options - options database, use NULL for default global database
2877: .  pre - string to prepend to each name or NULL
2878: .  name - the option one is seeking
2879: -  nmax - maximum number of values to retrieve

2881:    Output Parameters:
2882: +  dvalue - the scalar values to return
2883: .  nmax - actual number of values retreived
2884: -  set - PETSC_TRUE if found, else PETSC_FALSE

2886:    Level: beginner

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


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

2924: /*@C
2925:    PetscOptionsGetStringArray - Gets an array of string values for a particular
2926:    option in the database. The values must be separated with commas with
2927:    no intervening spaces.

2929:    Not Collective

2931:    Input Parameters:
2932: +  options - options database, use NULL for default global database
2933: .  pre - string to prepend to name or NULL
2934: .  name - the option one is seeking
2935: -  nmax - maximum number of strings

2937:    Output Parameters:
2938: +  strings - location to copy strings
2939: .  nmax - the number of strings found
2940: -  set - PETSC_TRUE if found, else PETSC_FALSE

2942:    Level: beginner

2944:    Notes:
2945:    The nmax parameter is used for both input and output.

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

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

2953: .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
2954:           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
2955:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
2956:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
2957:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
2958:           PetscOptionsFList(), PetscOptionsEList()
2959: @*/
2960: PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set)
2961: {
2962:   const char     *svalue;
2963:   char           *value;
2965:   PetscInt       n = 0;
2966:   PetscBool      flag;
2967:   PetscToken     token;


2974:   PetscOptionsFindPair(options,pre,name,&svalue,&flag);
2975:   if (!flag || !svalue)  { if (set) *set = PETSC_FALSE; *nmax = 0; return(0);}
2976:   if (set) *set = PETSC_TRUE;
2977:   PetscTokenCreate(svalue,',',&token);
2978:   PetscTokenFind(token,&value);
2979:   while (value && n < *nmax) {
2980:     PetscStrallocpy(value,&strings[n]);
2981:     PetscTokenFind(token,&value);
2982:     n++;
2983:   }
2984:   PetscTokenDestroy(&token);
2985:   *nmax = n;
2986:   return(0);
2987: }

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

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

2994:    Logically Collective

2996:    Input Parameters:
2997: +  pre - string to prepend to name or NULL
2998: .  oldname - the old, deprecated option
2999: .  newname - the new option, or NULL if option is purely removed
3000: .  version - a string describing the version of first deprecation, e.g. "3.9"
3001: -  info - additional information string, or NULL.

3003:    Options Database Keys:
3004: . -options_suppress_deprecated_warnings - do not print deprecation warnings

3006:    Notes:
3007:    Must be called between PetscOptionsBegin() (or PetscObjectOptionsBegin()) and PetscOptionsEnd().
3008:    Only the proces of rank zero that owns the PetscOptionsItems are argument (managed by PetscOptionsBegin() or
3009:    PetscObjectOptionsBegin() prints the information
3010:    If newname is provided, the old option is replaced. Otherwise, it remains
3011:    in the options database.
3012:    If an option is not replaced, the info argument should be used to advise the user
3013:    on how to proceed.
3014:    There is a limit on the length of the warning printed, so very long strings
3015:    provided as info may be truncated.

3017:    Level: developer

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

3021: @*/
3022: PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[])
3023: {
3024:   PetscErrorCode     ierr;
3025:   PetscBool          found,quiet;
3026:   const char         *value;
3027:   const char * const quietopt="-options_suppress_deprecated_warnings";
3028:   char               msg[4096];
3029:   char               *prefix = NULL;
3030:   PetscOptions       options = NULL;
3031:   MPI_Comm           comm = PETSC_COMM_SELF;

3036:   if (PetscOptionsObject) {
3037:     prefix  = PetscOptionsObject->prefix;
3038:     options = PetscOptionsObject->options;
3039:     comm    = PetscOptionsObject->comm;
3040:   }
3041:   PetscOptionsFindPair(options,prefix,oldname,&value,&found);
3042:   if (found) {
3043:     if (newname) {
3044:       if (prefix) {
3045:         PetscOptionsPrefixPush(options,prefix);
3046:       }
3047:       PetscOptionsSetValue(options,newname,value);
3048:       if (prefix) {
3049:         PetscOptionsPrefixPop(options);
3050:       }
3051:       PetscOptionsClearValue(options,oldname);
3052:     }
3053:     quiet = PETSC_FALSE;
3054:     PetscOptionsGetBool(options,NULL,quietopt,&quiet,NULL);
3055:     if (!quiet) {
3056:       PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");
3057:       PetscStrcat(msg,oldname);
3058:       PetscStrcat(msg," is deprecated as of version ");
3059:       PetscStrcat(msg,version);
3060:       PetscStrcat(msg," and will be removed in a future release.");
3061:       if (newname) {
3062:         PetscStrcat(msg," Please use the option ");
3063:         PetscStrcat(msg,newname);
3064:         PetscStrcat(msg," instead.");
3065:       }
3066:       if (info) {
3067:         PetscStrcat(msg," ");
3068:         PetscStrcat(msg,info);
3069:       }
3070:       PetscStrcat(msg," (Silence this warning with ");
3071:       PetscStrcat(msg,quietopt);
3072:       PetscStrcat(msg,")\n");
3073:       PetscPrintf(comm,msg);
3074:     }
3075:   }
3076:   return(0);
3077: }