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