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