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