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