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 .keywords: set, options, database 1613 @*/ 1614 PetscErrorCode PetscOptionsSetFromOptions(PetscOptions options) 1615 { 1616 PetscBool flgc = PETSC_FALSE,flgm; 1617 PetscErrorCode ierr; 1618 char monfilename[PETSC_MAX_PATH_LEN]; 1619 PetscViewer monviewer; 1620 1621 PetscFunctionBegin; 1622 /* 1623 The options argument is currently ignored since we currently maintain only a single options database 1624 1625 options = options ? options : defaultoptions; 1626 */ 1627 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for handling options","PetscOptions");CHKERRQ(ierr); 1628 ierr = PetscOptionsString("-options_monitor","Monitor options database","PetscOptionsMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flgm);CHKERRQ(ierr); 1629 ierr = PetscOptionsBool("-options_monitor_cancel","Cancel all options database monitors","PetscOptionsMonitorCancel",flgc,&flgc,NULL);CHKERRQ(ierr); 1630 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1631 if (flgm) { 1632 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,monfilename,&monviewer);CHKERRQ(ierr); 1633 ierr = PetscOptionsMonitorSet(PetscOptionsMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 1634 } 1635 if (flgc) { ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr); } 1636 PetscFunctionReturn(0); 1637 } 1638 1639 /*@C 1640 PetscOptionsMonitorDefault - Print all options set value events. 1641 1642 Logically Collective on PETSC_COMM_WORLD 1643 1644 Input Parameters: 1645 + name - option name string 1646 . value - option value string 1647 - ctx - an ASCII viewer 1648 1649 Level: intermediate 1650 1651 .keywords: PetscOptions, default, monitor 1652 1653 .seealso: PetscOptionsMonitorSet() 1654 @*/ 1655 PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx) 1656 { 1657 PetscErrorCode ierr; 1658 PetscViewer viewer = (PetscViewer)ctx; 1659 1660 PetscFunctionBegin; 1661 if (!value) { 1662 ierr = PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);CHKERRQ(ierr); 1663 } else if (!value[0]) { 1664 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);CHKERRQ(ierr); 1665 } else { 1666 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);CHKERRQ(ierr); 1667 } 1668 PetscFunctionReturn(0); 1669 } 1670 1671 /*@C 1672 PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that 1673 modified the PETSc options database. 1674 1675 Not Collective 1676 1677 Input Parameters: 1678 + monitor - pointer to function (if this is NULL, it turns off monitoring 1679 . mctx - [optional] context for private data for the 1680 monitor routine (use NULL if no context is desired) 1681 - monitordestroy - [optional] routine that frees monitor context 1682 (may be NULL) 1683 1684 Calling Sequence of monitor: 1685 $ monitor (const char name[], const char value[], void *mctx) 1686 1687 + name - option name string 1688 . value - option value string 1689 - mctx - optional monitoring context, as set by PetscOptionsMonitorSet() 1690 1691 Options Database Keys: 1692 + -options_monitor - sets PetscOptionsMonitorDefault() 1693 - -options_monitor_cancel - cancels all monitors that have 1694 been hardwired into a code by 1695 calls to PetscOptionsMonitorSet(), but 1696 does not cancel those set via 1697 the options database. 1698 1699 Notes: 1700 The default is to do nothing. To print the name and value of options 1701 being inserted into the database, use PetscOptionsMonitorDefault() as the monitoring routine, 1702 with a null monitoring context. 1703 1704 Several different monitoring routines may be set by calling 1705 PetscOptionsMonitorSet() multiple times; all will be called in the 1706 order in which they were set. 1707 1708 Level: beginner 1709 1710 .keywords: PetscOptions, set, monitor 1711 1712 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorCancel() 1713 @*/ 1714 PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 1715 { 1716 PetscOptions options = defaultoptions; 1717 1718 PetscFunctionBegin; 1719 if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set"); 1720 options->monitor[options->numbermonitors] = monitor; 1721 options->monitordestroy[options->numbermonitors] = monitordestroy; 1722 options->monitorcontext[options->numbermonitors++] = (void*)mctx; 1723 PetscFunctionReturn(0); 1724 } 1725 1726 /*@ 1727 PetscOptionsMonitorCancel - Clears all monitors for a PetscOptions object. 1728 1729 Not Collective 1730 1731 Options Database Key: 1732 . -options_monitor_cancel - Cancels all monitors that have 1733 been hardwired into a code by calls to PetscOptionsMonitorSet(), 1734 but does not cancel those set via the options database. 1735 1736 Level: intermediate 1737 1738 .keywords: PetscOptions, set, monitor 1739 1740 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorSet() 1741 @*/ 1742 PetscErrorCode PetscOptionsMonitorCancel(void) 1743 { 1744 PetscErrorCode ierr; 1745 PetscInt i; 1746 PetscOptions options = defaultoptions; 1747 1748 PetscFunctionBegin; 1749 for (i=0; i<options->numbermonitors; i++) { 1750 if (options->monitordestroy[i]) { 1751 ierr = (*options->monitordestroy[i])(&options->monitorcontext[i]);CHKERRQ(ierr); 1752 } 1753 } 1754 options->numbermonitors = 0; 1755 PetscFunctionReturn(0); 1756 } 1757 1758 /* 1759 PetscOptionsStringToBool - Converts string to PetscBool , handles cases like "yes", "no", "true", "false", "0", "1", "off", "on". 1760 */ 1761 PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a) 1762 { 1763 PetscBool istrue,isfalse; 1764 size_t len; 1765 PetscErrorCode ierr; 1766 1767 PetscFunctionBegin; 1768 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1769 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Character string of length zero has no logical value"); 1770 ierr = PetscStrcasecmp(value,"TRUE",&istrue);CHKERRQ(ierr); 1771 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1772 ierr = PetscStrcasecmp(value,"YES",&istrue);CHKERRQ(ierr); 1773 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1774 ierr = PetscStrcasecmp(value,"1",&istrue);CHKERRQ(ierr); 1775 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1776 ierr = PetscStrcasecmp(value,"on",&istrue);CHKERRQ(ierr); 1777 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1778 ierr = PetscStrcasecmp(value,"FALSE",&isfalse);CHKERRQ(ierr); 1779 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1780 ierr = PetscStrcasecmp(value,"NO",&isfalse);CHKERRQ(ierr); 1781 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1782 ierr = PetscStrcasecmp(value,"0",&isfalse);CHKERRQ(ierr); 1783 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1784 ierr = PetscStrcasecmp(value,"off",&isfalse);CHKERRQ(ierr); 1785 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1786 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value); 1787 } 1788 1789 /* 1790 PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide" 1791 */ 1792 PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a) 1793 { 1794 PetscErrorCode ierr; 1795 size_t len; 1796 PetscBool decide,tdefault,mouse; 1797 1798 PetscFunctionBegin; 1799 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1800 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 1801 1802 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);CHKERRQ(ierr); 1803 if (!tdefault) { 1804 ierr = PetscStrcasecmp(name,"DEFAULT",&tdefault);CHKERRQ(ierr); 1805 } 1806 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&decide);CHKERRQ(ierr); 1807 if (!decide) { 1808 ierr = PetscStrcasecmp(name,"DECIDE",&decide);CHKERRQ(ierr); 1809 } 1810 ierr = PetscStrcasecmp(name,"mouse",&mouse);CHKERRQ(ierr); 1811 1812 if (tdefault) *a = PETSC_DEFAULT; 1813 else if (decide) *a = PETSC_DECIDE; 1814 else if (mouse) *a = -1; 1815 else { 1816 char *endptr; 1817 long strtolval; 1818 1819 strtolval = strtol(name,&endptr,10); 1820 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); 1821 1822 #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL) 1823 (void) strtolval; 1824 *a = atoll(name); 1825 #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64) 1826 (void) strtolval; 1827 *a = _atoi64(name); 1828 #else 1829 *a = (PetscInt)strtolval; 1830 #endif 1831 } 1832 PetscFunctionReturn(0); 1833 } 1834 1835 #if defined(PETSC_USE_REAL___FLOAT128) 1836 #include <quadmath.h> 1837 #endif 1838 1839 static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr) 1840 { 1841 PetscFunctionBegin; 1842 #if defined(PETSC_USE_REAL___FLOAT128) 1843 *a = strtoflt128(name,endptr); 1844 #else 1845 *a = (PetscReal)strtod(name,endptr); 1846 #endif 1847 PetscFunctionReturn(0); 1848 } 1849 1850 static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary) 1851 { 1852 PetscBool hasi = PETSC_FALSE; 1853 char *ptr; 1854 PetscReal strtoval; 1855 PetscErrorCode ierr; 1856 1857 PetscFunctionBegin; 1858 ierr = PetscStrtod(name,&strtoval,&ptr);CHKERRQ(ierr); 1859 if (ptr == name) { 1860 strtoval = 1.; 1861 hasi = PETSC_TRUE; 1862 if (name[0] == 'i') { 1863 ptr++; 1864 } else if (name[0] == '+' && name[1] == 'i') { 1865 ptr += 2; 1866 } else if (name[0] == '-' && name[1] == 'i') { 1867 strtoval = -1.; 1868 ptr += 2; 1869 } 1870 } else if (*ptr == 'i') { 1871 hasi = PETSC_TRUE; 1872 ptr++; 1873 } 1874 *endptr = ptr; 1875 *isImaginary = hasi; 1876 if (hasi) { 1877 #if !defined(PETSC_USE_COMPLEX) 1878 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name); 1879 #else 1880 *a = PetscCMPLX(0.,strtoval); 1881 #endif 1882 } else { 1883 *a = strtoval; 1884 } 1885 PetscFunctionReturn(0); 1886 } 1887 1888 /* 1889 Converts a string to PetscReal value. Handles special cases like "default" and "decide" 1890 */ 1891 PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a) 1892 { 1893 size_t len; 1894 PetscBool match; 1895 char *endptr; 1896 PetscErrorCode ierr; 1897 1898 PetscFunctionBegin; 1899 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1900 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"String of length zero has no numerical value"); 1901 1902 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&match);CHKERRQ(ierr); 1903 if (!match) { 1904 ierr = PetscStrcasecmp(name,"DEFAULT",&match);CHKERRQ(ierr); 1905 } 1906 if (match) {*a = PETSC_DEFAULT; PetscFunctionReturn(0);} 1907 1908 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&match);CHKERRQ(ierr); 1909 if (!match) { 1910 ierr = PetscStrcasecmp(name,"DECIDE",&match);CHKERRQ(ierr); 1911 } 1912 if (match) {*a = PETSC_DECIDE; PetscFunctionReturn(0);} 1913 1914 ierr = PetscStrtod(name,a,&endptr);CHKERRQ(ierr); 1915 if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name); 1916 PetscFunctionReturn(0); 1917 } 1918 1919 PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a) 1920 { 1921 PetscBool imag1; 1922 size_t len; 1923 PetscScalar val = 0.; 1924 char *ptr = NULL; 1925 PetscErrorCode ierr; 1926 1927 PetscFunctionBegin; 1928 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1929 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 1930 ierr = PetscStrtoz(name,&val,&ptr,&imag1);CHKERRQ(ierr); 1931 #if defined(PETSC_USE_COMPLEX) 1932 if ((size_t) (ptr - name) < len) { 1933 PetscBool imag2; 1934 PetscScalar val2; 1935 1936 ierr = PetscStrtoz(ptr,&val2,&ptr,&imag2);CHKERRQ(ierr); 1937 if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name); 1938 val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2)); 1939 } 1940 #endif 1941 if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name); 1942 *a = val; 1943 PetscFunctionReturn(0); 1944 } 1945 1946 /*@C 1947 PetscOptionsGetBool - Gets the Logical (true or false) value for a particular 1948 option in the database. 1949 1950 Not Collective 1951 1952 Input Parameters: 1953 + options - options database, use NULL for default global database 1954 . pre - the string to prepend to the name or NULL 1955 - name - the option one is seeking 1956 1957 Output Parameter: 1958 + ivalue - the logical value to return 1959 - set - PETSC_TRUE if found, else PETSC_FALSE 1960 1961 Level: beginner 1962 1963 Notes: 1964 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 1965 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 1966 1967 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 1968 is equivalent to -requested_bool true 1969 1970 If the user does not supply the option at all ivalue is NOT changed. Thus 1971 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 1972 1973 Concepts: options database^has logical 1974 1975 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1976 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(), 1977 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1978 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1979 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1980 PetscOptionsFList(), PetscOptionsEList() 1981 @*/ 1982 PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set) 1983 { 1984 const char *value; 1985 PetscBool flag; 1986 PetscErrorCode ierr; 1987 1988 PetscFunctionBegin; 1989 PetscValidCharPointer(name,3); 1990 if (ivalue) PetscValidIntPointer(ivalue,4); 1991 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 1992 if (flag) { 1993 if (set) *set = PETSC_TRUE; 1994 if (!value) { 1995 if (ivalue) *ivalue = PETSC_TRUE; 1996 } else { 1997 ierr = PetscOptionsStringToBool(value, &flag);CHKERRQ(ierr); 1998 if (ivalue) *ivalue = flag; 1999 } 2000 } else { 2001 if (set) *set = PETSC_FALSE; 2002 } 2003 PetscFunctionReturn(0); 2004 } 2005 2006 /*@C 2007 PetscOptionsGetEList - Puts a list of option values that a single one may be selected from 2008 2009 Not Collective 2010 2011 Input Parameters: 2012 + options - options database, use NULL for default global database 2013 . pre - the string to prepend to the name or NULL 2014 . opt - option name 2015 . list - the possible choices (one of these must be selected, anything else is invalid) 2016 - ntext - number of choices 2017 2018 Output Parameter: 2019 + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed) 2020 - set - PETSC_TRUE if found, else PETSC_FALSE 2021 2022 Level: intermediate 2023 2024 Notes: 2025 If the user does not supply the option value is NOT changed. Thus 2026 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2027 2028 See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 2029 2030 Concepts: options database^list 2031 2032 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2033 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2034 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2035 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2036 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2037 PetscOptionsFList(), PetscOptionsEList() 2038 @*/ 2039 PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set) 2040 { 2041 PetscErrorCode ierr; 2042 size_t alen,len = 0; 2043 char *svalue; 2044 PetscBool aset,flg = PETSC_FALSE; 2045 PetscInt i; 2046 2047 PetscFunctionBegin; 2048 PetscValidCharPointer(opt,3); 2049 for (i=0; i<ntext; i++) { 2050 ierr = PetscStrlen(list[i],&alen);CHKERRQ(ierr); 2051 if (alen > len) len = alen; 2052 } 2053 len += 5; /* a little extra space for user mistypes */ 2054 ierr = PetscMalloc1(len,&svalue);CHKERRQ(ierr); 2055 ierr = PetscOptionsGetString(options,pre,opt,svalue,len,&aset);CHKERRQ(ierr); 2056 if (aset) { 2057 ierr = PetscEListFind(ntext,list,svalue,value,&flg);CHKERRQ(ierr); 2058 if (!flg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s",svalue,pre ? pre : "",opt+1); 2059 if (set) *set = PETSC_TRUE; 2060 } else if (set) *set = PETSC_FALSE; 2061 ierr = PetscFree(svalue);CHKERRQ(ierr); 2062 PetscFunctionReturn(0); 2063 } 2064 2065 /*@C 2066 PetscOptionsGetEnum - Gets the enum value for a particular option in the database. 2067 2068 Not Collective 2069 2070 Input Parameters: 2071 + options - options database, use NULL for default global database 2072 . pre - option prefix or NULL 2073 . opt - option name 2074 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2075 - defaultv - the default (current) value 2076 2077 Output Parameter: 2078 + value - the value to return 2079 - set - PETSC_TRUE if found, else PETSC_FALSE 2080 2081 Level: beginner 2082 2083 Concepts: options database 2084 2085 Notes: 2086 If the user does not supply the option value is NOT changed. Thus 2087 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2088 2089 List is usually something like PCASMTypes or some other predefined list of enum names 2090 2091 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2092 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2093 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2094 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2095 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2096 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2097 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2098 @*/ 2099 PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set) 2100 { 2101 PetscErrorCode ierr; 2102 PetscInt ntext = 0,tval; 2103 PetscBool fset; 2104 2105 PetscFunctionBegin; 2106 PetscValidCharPointer(opt,3); 2107 while (list[ntext++]) { 2108 if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 2109 } 2110 if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 2111 ntext -= 3; 2112 ierr = PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);CHKERRQ(ierr); 2113 /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 2114 if (fset) *value = (PetscEnum)tval; 2115 if (set) *set = fset; 2116 PetscFunctionReturn(0); 2117 } 2118 2119 /*@C 2120 PetscOptionsGetInt - Gets the integer value for a particular option in the database. 2121 2122 Not Collective 2123 2124 Input Parameters: 2125 + options - options database, use NULL for default global database 2126 . pre - the string to prepend to the name or NULL 2127 - name - the option one is seeking 2128 2129 Output Parameter: 2130 + ivalue - the integer value to return 2131 - set - PETSC_TRUE if found, else PETSC_FALSE 2132 2133 Level: beginner 2134 2135 Notes: 2136 If the user does not supply the option ivalue is NOT changed. Thus 2137 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2138 2139 Concepts: options database^has int 2140 2141 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 2142 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2143 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2144 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2145 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2146 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2147 PetscOptionsFList(), PetscOptionsEList() 2148 @*/ 2149 PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set) 2150 { 2151 const char *value; 2152 PetscErrorCode ierr; 2153 PetscBool flag; 2154 2155 PetscFunctionBegin; 2156 PetscValidCharPointer(name,3); 2157 PetscValidIntPointer(ivalue,4); 2158 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2159 if (flag) { 2160 if (!value) { 2161 if (set) *set = PETSC_FALSE; 2162 } else { 2163 if (set) *set = PETSC_TRUE; 2164 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2165 } 2166 } else { 2167 if (set) *set = PETSC_FALSE; 2168 } 2169 PetscFunctionReturn(0); 2170 } 2171 2172 /*@C 2173 PetscOptionsGetReal - Gets the double precision value for a particular 2174 option in the database. 2175 2176 Not Collective 2177 2178 Input Parameters: 2179 + options - options database, use NULL for default global database 2180 . pre - string to prepend to each name or NULL 2181 - name - the option one is seeking 2182 2183 Output Parameter: 2184 + dvalue - the double value to return 2185 - set - PETSC_TRUE if found, PETSC_FALSE if not found 2186 2187 Notes: 2188 If the user does not supply the option dvalue is NOT changed. Thus 2189 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2190 2191 Level: beginner 2192 2193 Concepts: options database^has double 2194 2195 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2196 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(), 2197 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2198 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2199 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2200 PetscOptionsFList(), PetscOptionsEList() 2201 @*/ 2202 PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set) 2203 { 2204 const char *value; 2205 PetscBool flag; 2206 PetscErrorCode ierr; 2207 2208 PetscFunctionBegin; 2209 PetscValidCharPointer(name,3); 2210 PetscValidRealPointer(dvalue,4); 2211 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2212 if (flag) { 2213 if (!value) { 2214 if (set) *set = PETSC_FALSE; 2215 } else { 2216 if (set) *set = PETSC_TRUE; 2217 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2218 } 2219 } else { 2220 if (set) *set = PETSC_FALSE; 2221 } 2222 PetscFunctionReturn(0); 2223 } 2224 2225 /*@C 2226 PetscOptionsGetScalar - Gets the scalar value for a particular 2227 option in the database. 2228 2229 Not Collective 2230 2231 Input Parameters: 2232 + options - options database, use NULL for default global database 2233 . pre - string to prepend to each name or NULL 2234 - name - the option one is seeking 2235 2236 Output Parameter: 2237 + dvalue - the double value to return 2238 - set - PETSC_TRUE if found, else PETSC_FALSE 2239 2240 Level: beginner 2241 2242 Usage: 2243 A complex number 2+3i must be specified with NO spaces 2244 2245 Notes: 2246 If the user does not supply the option dvalue is NOT changed. Thus 2247 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2248 2249 Concepts: options database^has scalar 2250 2251 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2252 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2253 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2254 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2255 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2256 PetscOptionsFList(), PetscOptionsEList() 2257 @*/ 2258 PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set) 2259 { 2260 const char *value; 2261 PetscBool flag; 2262 PetscErrorCode ierr; 2263 2264 PetscFunctionBegin; 2265 PetscValidCharPointer(name,3); 2266 PetscValidScalarPointer(dvalue,4); 2267 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2268 if (flag) { 2269 if (!value) { 2270 if (set) *set = PETSC_FALSE; 2271 } else { 2272 #if !defined(PETSC_USE_COMPLEX) 2273 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2274 #else 2275 ierr = PetscOptionsStringToScalar(value,dvalue);CHKERRQ(ierr); 2276 #endif 2277 if (set) *set = PETSC_TRUE; 2278 } 2279 } else { /* flag */ 2280 if (set) *set = PETSC_FALSE; 2281 } 2282 PetscFunctionReturn(0); 2283 } 2284 2285 /*@C 2286 PetscOptionsGetString - Gets the string value for a particular option in 2287 the database. 2288 2289 Not Collective 2290 2291 Input Parameters: 2292 + options - options database, use NULL for default global database 2293 . pre - string to prepend to name or NULL 2294 . name - the option one is seeking 2295 - len - maximum length of the string including null termination 2296 2297 Output Parameters: 2298 + string - location to copy string 2299 - set - PETSC_TRUE if found, else PETSC_FALSE 2300 2301 Level: beginner 2302 2303 Fortran Note: 2304 The Fortran interface is slightly different from the C/C++ 2305 interface (len is not used). Sample usage in Fortran follows 2306 .vb 2307 character *20 string 2308 PetscErrorCode ierr 2309 PetscBool set 2310 call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr) 2311 .ve 2312 2313 Notes: 2314 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 2315 2316 If the user does not use the option then the string is not changed. Thus 2317 you should ALWAYS initialize the string if you access it without first checking if the set flag is true. 2318 2319 Concepts: options database^string 2320 2321 Note: 2322 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). 2323 2324 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2325 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2326 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2327 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2328 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2329 PetscOptionsFList(), PetscOptionsEList() 2330 @*/ 2331 PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set) 2332 { 2333 const char *value; 2334 PetscBool flag; 2335 PetscErrorCode ierr; 2336 2337 PetscFunctionBegin; 2338 PetscValidCharPointer(name,3); 2339 PetscValidCharPointer(string,4); 2340 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2341 if (!flag) { 2342 if (set) *set = PETSC_FALSE; 2343 } else { 2344 if (set) *set = PETSC_TRUE; 2345 if (value) { 2346 ierr = PetscStrncpy(string,value,len);CHKERRQ(ierr); 2347 } else { 2348 ierr = PetscMemzero(string,len);CHKERRQ(ierr); 2349 } 2350 } 2351 PetscFunctionReturn(0); 2352 } 2353 2354 char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[]) 2355 { 2356 const char *value; 2357 PetscBool flag; 2358 PetscErrorCode ierr; 2359 2360 PetscFunctionBegin; 2361 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(0); 2362 if (flag) PetscFunctionReturn((char*)value); 2363 else PetscFunctionReturn(0); 2364 } 2365 2366 /*@C 2367 PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular 2368 option in the database. The values must be separated with commas with 2369 no intervening spaces. 2370 2371 Not Collective 2372 2373 Input Parameters: 2374 + options - options database, use NULL for default global database 2375 . pre - string to prepend to each name or NULL 2376 . name - the option one is seeking 2377 - nmax - maximum number of values to retrieve 2378 2379 Output Parameter: 2380 + dvalue - the integer values to return 2381 . nmax - actual number of values retreived 2382 - set - PETSC_TRUE if found, else PETSC_FALSE 2383 2384 Level: beginner 2385 2386 Concepts: options database^array of ints 2387 2388 Notes: 2389 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 2390 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 2391 2392 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2393 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2394 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2395 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2396 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2397 PetscOptionsFList(), PetscOptionsEList() 2398 @*/ 2399 PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set) 2400 { 2401 const char *svalue; 2402 char *value; 2403 PetscErrorCode ierr; 2404 PetscInt n = 0; 2405 PetscBool flag; 2406 PetscToken token; 2407 2408 PetscFunctionBegin; 2409 PetscValidCharPointer(name,3); 2410 PetscValidIntPointer(dvalue,4); 2411 PetscValidIntPointer(nmax,5); 2412 2413 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2414 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2415 if (set) *set = PETSC_TRUE; 2416 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2417 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2418 while (value && n < *nmax) { 2419 ierr = PetscOptionsStringToBool(value,dvalue);CHKERRQ(ierr); 2420 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2421 dvalue++; 2422 n++; 2423 } 2424 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2425 *nmax = n; 2426 PetscFunctionReturn(0); 2427 } 2428 2429 /*@C 2430 PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database. 2431 2432 Not Collective 2433 2434 Input Parameters: 2435 + options - options database, use NULL for default global database 2436 . pre - option prefix or NULL 2437 . name - option name 2438 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2439 - nmax - maximum number of values to retrieve 2440 2441 Output Parameters: 2442 + ivalue - the enum values to return 2443 . nmax - actual number of values retreived 2444 - set - PETSC_TRUE if found, else PETSC_FALSE 2445 2446 Level: beginner 2447 2448 Concepts: options database 2449 2450 Notes: 2451 The array must be passed as a comma separated list. 2452 2453 There must be no intervening spaces between the values. 2454 2455 list is usually something like PCASMTypes or some other predefined list of enum names. 2456 2457 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2458 PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2459 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(), 2460 PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(), 2461 PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2462 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2463 @*/ 2464 PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set) 2465 { 2466 const char *svalue; 2467 char *value; 2468 PetscInt n = 0; 2469 PetscEnum evalue; 2470 PetscBool flag; 2471 PetscToken token; 2472 PetscErrorCode ierr; 2473 2474 PetscFunctionBegin; 2475 PetscValidCharPointer(name,3); 2476 PetscValidPointer(list,4); 2477 PetscValidPointer(ivalue,5); 2478 PetscValidIntPointer(nmax,6); 2479 2480 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2481 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2482 if (set) *set = PETSC_TRUE; 2483 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2484 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2485 while (value && n < *nmax) { 2486 ierr = PetscEnumFind(list,value,&evalue,&flag);CHKERRQ(ierr); 2487 if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1); 2488 ivalue[n++] = evalue; 2489 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2490 } 2491 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2492 *nmax = n; 2493 PetscFunctionReturn(0); 2494 } 2495 2496 /*@C 2497 PetscOptionsGetIntArray - Gets an array of integer values for a particular 2498 option in the database. 2499 2500 Not Collective 2501 2502 Input Parameters: 2503 + options - options database, use NULL for default global database 2504 . pre - string to prepend to each name or NULL 2505 . name - the option one is seeking 2506 - nmax - maximum number of values to retrieve 2507 2508 Output Parameter: 2509 + ivalue - the integer values to return 2510 . nmax - actual number of values retreived 2511 - set - PETSC_TRUE if found, else PETSC_FALSE 2512 2513 Level: beginner 2514 2515 Notes: 2516 The array can be passed as 2517 a comma separated list: 0,1,2,3,4,5,6,7 2518 a range (start-end+1): 0-8 2519 a range with given increment (start-end+1:inc): 0-7:2 2520 a combination of values and ranges separated by commas: 0,1-8,8-15:2 2521 2522 There must be no intervening spaces between the values. 2523 2524 Concepts: options database^array of ints 2525 2526 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2527 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2528 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2529 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2530 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2531 PetscOptionsFList(), PetscOptionsEList() 2532 @*/ 2533 PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set) 2534 { 2535 const char *svalue; 2536 char *value; 2537 PetscErrorCode ierr; 2538 PetscInt n = 0,i,j,start,end,inc,nvalues; 2539 size_t len; 2540 PetscBool flag,foundrange; 2541 PetscToken token; 2542 2543 PetscFunctionBegin; 2544 PetscValidCharPointer(name,3); 2545 PetscValidIntPointer(ivalue,4); 2546 PetscValidIntPointer(nmax,5); 2547 2548 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2549 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2550 if (set) *set = PETSC_TRUE; 2551 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2552 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2553 while (value && n < *nmax) { 2554 /* look for form d-D where d and D are integers */ 2555 foundrange = PETSC_FALSE; 2556 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 2557 if (value[0] == '-') i=2; 2558 else i=1; 2559 for (;i<(int)len; i++) { 2560 if (value[i] == '-') { 2561 if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 2562 value[i] = 0; 2563 2564 ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 2565 inc = 1; 2566 j = i+1; 2567 for (;j<(int)len; j++) { 2568 if (value[j] == ':') { 2569 value[j] = 0; 2570 2571 ierr = PetscOptionsStringToInt(value+j+1,&inc);CHKERRQ(ierr); 2572 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); 2573 break; 2574 } 2575 } 2576 ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 2577 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); 2578 nvalues = (end-start)/inc + (end-start)%inc; 2579 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); 2580 for (;start<end; start+=inc) { 2581 *ivalue = start; ivalue++;n++; 2582 } 2583 foundrange = PETSC_TRUE; 2584 break; 2585 } 2586 } 2587 if (!foundrange) { 2588 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2589 ivalue++; 2590 n++; 2591 } 2592 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2593 } 2594 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2595 *nmax = n; 2596 PetscFunctionReturn(0); 2597 } 2598 2599 /*@C 2600 PetscOptionsGetRealArray - Gets an array of double precision values for a 2601 particular option in the database. The values must be separated with 2602 commas with no intervening spaces. 2603 2604 Not Collective 2605 2606 Input Parameters: 2607 + options - options database, use NULL for default global database 2608 . pre - string to prepend to each name or NULL 2609 . name - the option one is seeking 2610 - nmax - maximum number of values to retrieve 2611 2612 Output Parameters: 2613 + dvalue - the double values to return 2614 . nmax - actual number of values retreived 2615 - set - PETSC_TRUE if found, else PETSC_FALSE 2616 2617 Level: beginner 2618 2619 Concepts: options database^array of doubles 2620 2621 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2622 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2623 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2624 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2625 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2626 PetscOptionsFList(), PetscOptionsEList() 2627 @*/ 2628 PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set) 2629 { 2630 const char *svalue; 2631 char *value; 2632 PetscErrorCode ierr; 2633 PetscInt n = 0; 2634 PetscBool flag; 2635 PetscToken token; 2636 2637 PetscFunctionBegin; 2638 PetscValidCharPointer(name,3); 2639 PetscValidRealPointer(dvalue,4); 2640 PetscValidIntPointer(nmax,5); 2641 2642 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2643 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2644 if (set) *set = PETSC_TRUE; 2645 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2646 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2647 while (value && n < *nmax) { 2648 ierr = PetscOptionsStringToReal(value,dvalue++);CHKERRQ(ierr); 2649 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2650 n++; 2651 } 2652 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2653 *nmax = n; 2654 PetscFunctionReturn(0); 2655 } 2656 2657 /*@C 2658 PetscOptionsGetScalarArray - Gets an array of scalars for a 2659 particular option in the database. The values must be separated with 2660 commas with no intervening spaces. 2661 2662 Not Collective 2663 2664 Input Parameters: 2665 + options - options database, use NULL for default global database 2666 . pre - string to prepend to each name or NULL 2667 . name - the option one is seeking 2668 - nmax - maximum number of values to retrieve 2669 2670 Output Parameters: 2671 + dvalue - the scalar values to return 2672 . nmax - actual number of values retreived 2673 - set - PETSC_TRUE if found, else PETSC_FALSE 2674 2675 Level: beginner 2676 2677 Concepts: options database^array of doubles 2678 2679 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2680 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2681 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2682 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2683 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2684 PetscOptionsFList(), PetscOptionsEList() 2685 @*/ 2686 PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],PetscInt *nmax,PetscBool *set) 2687 { 2688 const char *svalue; 2689 char *value; 2690 PetscErrorCode ierr; 2691 PetscInt n = 0; 2692 PetscBool flag; 2693 PetscToken token; 2694 2695 PetscFunctionBegin; 2696 PetscValidCharPointer(name,3); 2697 PetscValidRealPointer(dvalue,4); 2698 PetscValidIntPointer(nmax,5); 2699 2700 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2701 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2702 if (set) *set = PETSC_TRUE; 2703 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2704 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2705 while (value && n < *nmax) { 2706 ierr = PetscOptionsStringToScalar(value,dvalue++);CHKERRQ(ierr); 2707 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2708 n++; 2709 } 2710 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2711 *nmax = n; 2712 PetscFunctionReturn(0); 2713 } 2714 2715 /*@C 2716 PetscOptionsGetStringArray - Gets an array of string values for a particular 2717 option in the database. The values must be separated with commas with 2718 no intervening spaces. 2719 2720 Not Collective 2721 2722 Input Parameters: 2723 + options - options database, use NULL for default global database 2724 . pre - string to prepend to name or NULL 2725 . name - the option one is seeking 2726 - nmax - maximum number of strings 2727 2728 Output Parameter: 2729 + strings - location to copy strings 2730 - set - PETSC_TRUE if found, else PETSC_FALSE 2731 2732 Level: beginner 2733 2734 Notes: 2735 The user should pass in an array of pointers to char, to hold all the 2736 strings returned by this function. 2737 2738 The user is responsible for deallocating the strings that are 2739 returned. The Fortran interface for this routine is not supported. 2740 2741 Contributed by Matthew Knepley. 2742 2743 Concepts: options database^array of strings 2744 2745 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2746 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2747 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2748 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2749 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2750 PetscOptionsFList(), PetscOptionsEList() 2751 @*/ 2752 PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set) 2753 { 2754 const char *svalue; 2755 char *value; 2756 PetscErrorCode ierr; 2757 PetscInt n = 0; 2758 PetscBool flag; 2759 PetscToken token; 2760 2761 PetscFunctionBegin; 2762 PetscValidCharPointer(name,3); 2763 PetscValidPointer(strings,4); 2764 PetscValidIntPointer(nmax,5); 2765 2766 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2767 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2768 if (set) *set = PETSC_TRUE; 2769 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2770 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2771 while (value && n < *nmax) { 2772 ierr = PetscStrallocpy(value,&strings[n]);CHKERRQ(ierr); 2773 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2774 n++; 2775 } 2776 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2777 *nmax = n; 2778 PetscFunctionReturn(0); 2779 } 2780 2781 /*@C 2782 PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one 2783 2784 Prints a deprecation warning, unless an option is supplied to suppress. 2785 2786 Not Collective 2787 2788 Input Parameters: 2789 + pre - string to prepend to name or NULL 2790 . oldname - the old, deprecated option 2791 . newname - the new option, or NULL if option is purely removed 2792 . version - a string describing the version of first deprecation, e.g. "3.9" 2793 - info - additional information string, or NULL. 2794 2795 Options Database Keys: 2796 . -options_suppress_deprecated_warnings - do not print deprecation warnings 2797 2798 Notes: 2799 Must be called between PetscOptionsBegin() and PetscOptionsEnd(). 2800 If newname is provided, the old option is replaced. Otherwise, it remains 2801 in the options database. 2802 If an option is not replaced, the info argument should be used to advise the user 2803 on how to proceed. 2804 There is a limit on the length of the warning printed, so very long strings 2805 provided as info may be truncated. 2806 2807 Level: developer 2808 2809 .seealso: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsScalar(), PetscOptionsBool(), PetscOptionsString(), PetscOptionsSetValue() 2810 2811 @*/ 2812 PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[]) 2813 { 2814 PetscErrorCode ierr; 2815 PetscBool found,quiet; 2816 const char *value; 2817 const char * const quietopt="-options_suppress_deprecated_warnings"; 2818 char msg[4096]; 2819 2820 PetscFunctionBegin; 2821 PetscValidCharPointer(oldname,2); 2822 PetscValidCharPointer(version,4); 2823 2824 ierr = PetscOptionsFindPair(PetscOptionsObject->options,PetscOptionsObject->prefix,oldname,&value,&found);CHKERRQ(ierr); 2825 if (found) { 2826 if (newname) { 2827 if (PetscOptionsObject->prefix) { 2828 ierr = PetscOptionsPrefixPush(PetscOptionsObject->options,PetscOptionsObject->prefix);CHKERRQ(ierr); 2829 } 2830 ierr = PetscOptionsSetValue(PetscOptionsObject->options,newname,value);CHKERRQ(ierr); 2831 if (PetscOptionsObject->prefix) { 2832 ierr = PetscOptionsPrefixPop(PetscOptionsObject->options);CHKERRQ(ierr); 2833 } 2834 ierr = PetscOptionsClearValue(PetscOptionsObject->options,oldname);CHKERRQ(ierr); 2835 } 2836 quiet = PETSC_FALSE; 2837 ierr = PetscOptionsGetBool(PetscOptionsObject->options,NULL,quietopt,&quiet,NULL);CHKERRQ(ierr); 2838 if (!quiet) { 2839 ierr = PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");CHKERRQ(ierr); 2840 ierr = PetscStrcat(msg,oldname);CHKERRQ(ierr); 2841 ierr = PetscStrcat(msg," is deprecated as of version ");CHKERRQ(ierr); 2842 ierr = PetscStrcat(msg,version);CHKERRQ(ierr); 2843 ierr = PetscStrcat(msg," and will be removed in a future release.");CHKERRQ(ierr); 2844 if (newname) { 2845 ierr = PetscStrcat(msg," Please use the option ");CHKERRQ(ierr); 2846 ierr = PetscStrcat(msg,newname);CHKERRQ(ierr); 2847 ierr = PetscStrcat(msg," instead.");CHKERRQ(ierr); 2848 } 2849 if (info) { 2850 ierr = PetscStrcat(msg," ");CHKERRQ(ierr); 2851 ierr = PetscStrcat(msg,info);CHKERRQ(ierr); 2852 } 2853 ierr = PetscStrcat(msg," (Silence this warning with ");CHKERRQ(ierr); 2854 ierr = PetscStrcat(msg,quietopt);CHKERRQ(ierr); 2855 ierr = PetscStrcat(msg,")\n");CHKERRQ(ierr); 2856 ierr = PetscPrintf(PetscOptionsObject->comm,msg);CHKERRQ(ierr); 2857 } 2858 } 2859 PetscFunctionReturn(0); 2860 } 2861