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