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