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