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