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