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