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 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[],const char *value[],PetscBool *set) 1170 { 1171 char buf[MAXOPTNAME]; 1172 PetscErrorCode ierr; 1173 1174 PetscFunctionBegin; 1175 options = options ? options : defaultoptions; 1176 if (pre && pre[0] == '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre); 1177 if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name); 1178 1179 name++; /* skip starting dash */ 1180 1181 /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */ 1182 if (pre && pre[0]) { 1183 char *ptr = buf; 1184 if (name[0] == '-') { *ptr++ = '-'; name++; } 1185 ierr = PetscStrncpy(ptr,pre,sizeof(buf)+(size_t)(ptr-buf));CHKERRQ(ierr); 1186 ierr = PetscStrlcat(buf,name,sizeof(buf));CHKERRQ(ierr); 1187 name = buf; 1188 } 1189 1190 #if defined(PETSC_USE_DEBUG) 1191 { 1192 PetscBool valid; 1193 char key[MAXOPTNAME+1] = "-"; 1194 ierr = PetscStrncpy(key+1,name,sizeof(key)-1);CHKERRQ(ierr); 1195 ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr); 1196 if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name); 1197 } 1198 #endif 1199 1200 { /* slow search */ 1201 int i; 1202 size_t len; 1203 PetscBool match; 1204 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1205 for (i=0; i<options->N; i++) { 1206 ierr = PetscStrncmp(options->names[i],name,len,&match);CHKERRQ(ierr); 1207 if (match) { 1208 options->used[i] = PETSC_TRUE; 1209 if (value) *value = options->values[i]; 1210 if (set) *set = PETSC_TRUE; 1211 PetscFunctionReturn(0); 1212 } 1213 } 1214 } 1215 1216 if (set) *set = PETSC_FALSE; 1217 PetscFunctionReturn(0); 1218 } 1219 1220 /*@C 1221 PetscOptionsReject - Generates an error if a certain option is given. 1222 1223 Not Collective, but setting values on certain processors could cause problems 1224 for parallel objects looking for options. 1225 1226 Input Parameters: 1227 + options - options database, use NULL for default global database 1228 . pre - the option prefix (may be NULL) 1229 . name - the option name one is seeking 1230 - mess - error message (may be NULL) 1231 1232 Level: advanced 1233 1234 Concepts: options database^rejecting option 1235 1236 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(), 1237 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1238 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1239 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1240 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1241 PetscOptionsFList(), PetscOptionsEList() 1242 @*/ 1243 PetscErrorCode PetscOptionsReject(PetscOptions options,const char pre[],const char name[],const char mess[]) 1244 { 1245 PetscErrorCode ierr; 1246 PetscBool flag = PETSC_FALSE; 1247 1248 PetscFunctionBegin; 1249 ierr = PetscOptionsHasName(options,pre,name,&flag);CHKERRQ(ierr); 1250 if (flag) { 1251 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); 1252 else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s",pre?pre:"",name+1); 1253 } 1254 PetscFunctionReturn(0); 1255 } 1256 1257 /*@C 1258 PetscOptionsHasHelp - Determines whether the "-help" option is in the database. 1259 1260 Not Collective 1261 1262 Input Parameters: 1263 . options - options database, use NULL for default global database 1264 1265 Output Parameters: 1266 . set - PETSC_TRUE if found else PETSC_FALSE. 1267 1268 Level: advanced 1269 1270 Concepts: options database^help 1271 1272 .seealso: PetscOptionsHasName() 1273 @*/ 1274 PetscErrorCode PetscOptionsHasHelp(PetscOptions options,PetscBool *set) 1275 { 1276 PetscFunctionBegin; 1277 PetscValidPointer(set,2); 1278 options = options ? options : defaultoptions; 1279 *set = options->help; 1280 PetscFunctionReturn(0); 1281 } 1282 1283 /*@C 1284 PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or boolean, even 1285 its value is set to false. 1286 1287 Not Collective 1288 1289 Input Parameters: 1290 + options - options database, use NULL for default global database 1291 . pre - string to prepend to the name or NULL 1292 - name - the option one is seeking 1293 1294 Output Parameters: 1295 . set - PETSC_TRUE if found else PETSC_FALSE. 1296 1297 Level: beginner 1298 1299 Concepts: options database^has option name 1300 1301 Notes: 1302 Name cannot be simply "-h". 1303 1304 In many cases you probably want to use PetscOptionsGetBool() instead of calling this, to allowing toggling values. 1305 1306 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1307 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1308 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1309 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1310 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1311 PetscOptionsFList(), PetscOptionsEList() 1312 @*/ 1313 PetscErrorCode PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool *set) 1314 { 1315 const char *value; 1316 PetscErrorCode ierr; 1317 PetscBool flag; 1318 1319 PetscFunctionBegin; 1320 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 1321 if (set) *set = flag; 1322 PetscFunctionReturn(0); 1323 } 1324 1325 /*@C 1326 PetscOptionsGetAll - Lists all the options the program was run with in a single string. 1327 1328 Not Collective 1329 1330 Input Paramter: 1331 . options - the options database, use NULL for the default global database 1332 1333 Output Parameter: 1334 . copts - pointer where string pointer is stored 1335 1336 Notes: 1337 the array and each entry in the array should be freed with PetscFree() 1338 1339 Level: advanced 1340 1341 Concepts: options database^listing 1342 1343 .seealso: PetscOptionsAllUsed(), PetscOptionsView() 1344 @*/ 1345 PetscErrorCode PetscOptionsGetAll(PetscOptions options,char *copts[]) 1346 { 1347 PetscErrorCode ierr; 1348 PetscInt i; 1349 size_t len = 1,lent = 0; 1350 char *coptions = NULL; 1351 1352 PetscFunctionBegin; 1353 PetscValidPointer(copts,2); 1354 options = options ? options : defaultoptions; 1355 /* count the length of the required string */ 1356 for (i=0; i<options->N; i++) { 1357 ierr = PetscStrlen(options->names[i],&lent);CHKERRQ(ierr); 1358 len += 2 + lent; 1359 if (options->values[i]) { 1360 ierr = PetscStrlen(options->values[i],&lent);CHKERRQ(ierr); 1361 len += 1 + lent; 1362 } 1363 } 1364 ierr = PetscMalloc1(len,&coptions);CHKERRQ(ierr); 1365 coptions[0] = 0; 1366 for (i=0; i<options->N; i++) { 1367 ierr = PetscStrcat(coptions,"-");CHKERRQ(ierr); 1368 ierr = PetscStrcat(coptions,options->names[i]);CHKERRQ(ierr); 1369 ierr = PetscStrcat(coptions," ");CHKERRQ(ierr); 1370 if (options->values[i]) { 1371 ierr = PetscStrcat(coptions,options->values[i]);CHKERRQ(ierr); 1372 ierr = PetscStrcat(coptions," ");CHKERRQ(ierr); 1373 } 1374 } 1375 *copts = coptions; 1376 PetscFunctionReturn(0); 1377 } 1378 1379 /*@C 1380 PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database 1381 1382 Not Collective 1383 1384 Input Parameter: 1385 + options - options database, use NULL for default global database 1386 - name - string name of option 1387 1388 Output Parameter: 1389 . used - PETSC_TRUE if the option was used, otherwise false, including if option was not found in options database 1390 1391 Level: advanced 1392 1393 .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed() 1394 @*/ 1395 PetscErrorCode PetscOptionsUsed(PetscOptions options,const char *name,PetscBool *used) 1396 { 1397 PetscInt i; 1398 PetscErrorCode ierr; 1399 1400 PetscFunctionBegin; 1401 PetscValidCharPointer(name,2); 1402 PetscValidPointer(used,3); 1403 options = options ? options : defaultoptions; 1404 *used = PETSC_FALSE; 1405 for (i=0; i<options->N; i++) { 1406 ierr = PetscStrcmp(options->names[i],name,used);CHKERRQ(ierr); 1407 if (*used) { 1408 *used = options->used[i]; 1409 break; 1410 } 1411 } 1412 PetscFunctionReturn(0); 1413 } 1414 1415 /*@ 1416 PetscOptionsAllUsed - Returns a count of the number of options in the 1417 database that have never been selected. 1418 1419 Not Collective 1420 1421 Input Parameter: 1422 . options - options database, use NULL for default global database 1423 1424 Output Parameter: 1425 . N - count of options not used 1426 1427 Level: advanced 1428 1429 .seealso: PetscOptionsView() 1430 @*/ 1431 PetscErrorCode PetscOptionsAllUsed(PetscOptions options,PetscInt *N) 1432 { 1433 PetscInt i,n = 0; 1434 1435 PetscFunctionBegin; 1436 PetscValidIntPointer(N,2); 1437 options = options ? options : defaultoptions; 1438 for (i=0; i<options->N; i++) { 1439 if (!options->used[i]) n++; 1440 } 1441 *N = n; 1442 PetscFunctionReturn(0); 1443 } 1444 1445 /*@ 1446 PetscOptionsLeft - Prints to screen any options that were set and never used. 1447 1448 Not Collective 1449 1450 Input Parameter: 1451 . options - options database; use NULL for default global database 1452 1453 Options Database Key: 1454 . -options_left - Activates OptionsAllUsed() within PetscFinalize() 1455 1456 Level: advanced 1457 1458 .seealso: PetscOptionsAllUsed() 1459 @*/ 1460 PetscErrorCode PetscOptionsLeft(PetscOptions options) 1461 { 1462 PetscErrorCode ierr; 1463 PetscInt i; 1464 1465 PetscFunctionBegin; 1466 options = options ? options : defaultoptions; 1467 for (i=0; i<options->N; i++) { 1468 if (!options->used[i]) { 1469 if (options->values[i]) { 1470 ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",options->names[i],options->values[i]);CHKERRQ(ierr); 1471 } else { 1472 ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",options->names[i]);CHKERRQ(ierr); 1473 } 1474 } 1475 } 1476 PetscFunctionReturn(0); 1477 } 1478 1479 /*@C 1480 PetscOptionsLeftGet - Returns all options that were set and never used. 1481 1482 Not Collective 1483 1484 Input Parameter: 1485 . options - options database, use NULL for default global database 1486 1487 Output Parameter: 1488 . N - count of options not used 1489 . names - names of options not used 1490 . values - values of options not used 1491 1492 Level: advanced 1493 1494 Notes: 1495 Users should call PetscOptionsLeftRestore() to free the memory allocated in this routine 1496 1497 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft() 1498 @*/ 1499 PetscErrorCode PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[]) 1500 { 1501 PetscErrorCode ierr; 1502 PetscInt i,n; 1503 1504 PetscFunctionBegin; 1505 if (N) PetscValidIntPointer(N,2); 1506 if (names) PetscValidPointer(names,3); 1507 if (values) PetscValidPointer(values,4); 1508 options = options ? options : defaultoptions; 1509 1510 /* The number of unused PETSc options */ 1511 n = 0; 1512 for (i=0; i<options->N; i++) { 1513 if (!options->used[i]) n++; 1514 } 1515 if (N) { *N = n; } 1516 if (names) { ierr = PetscMalloc1(n,names);CHKERRQ(ierr); } 1517 if (values) { ierr = PetscMalloc1(n,values);CHKERRQ(ierr); } 1518 1519 n = 0; 1520 if (names || values) { 1521 for (i=0; i<options->N; i++) { 1522 if (!options->used[i]) { 1523 if (names) (*names)[n] = options->names[i]; 1524 if (values) (*values)[n] = options->values[i]; 1525 n++; 1526 } 1527 } 1528 } 1529 PetscFunctionReturn(0); 1530 } 1531 1532 /*@C 1533 PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using PetscOptionsLeftGet. 1534 1535 Not Collective 1536 1537 Input Parameter: 1538 . options - options database, use NULL for default global database 1539 . names - names of options not used 1540 . values - values of options not used 1541 1542 Level: advanced 1543 1544 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet() 1545 @*/ 1546 PetscErrorCode PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[]) 1547 { 1548 PetscErrorCode ierr; 1549 1550 PetscFunctionBegin; 1551 if (N) PetscValidIntPointer(N,2); 1552 if (names) PetscValidPointer(names,3); 1553 if (values) PetscValidPointer(values,4); 1554 if (N) { *N = 0; } 1555 if (names) { ierr = PetscFree(*names);CHKERRQ(ierr); } 1556 if (values) { ierr = PetscFree(*values);CHKERRQ(ierr); } 1557 PetscFunctionReturn(0); 1558 } 1559 1560 /*@C 1561 PetscOptionsSetFromOptions - Sets options related to the handling of options in PETSc 1562 1563 Collective on PETSC_COMM_WORLD 1564 1565 Input Parameter: 1566 . options - options database, use NULL for default global database 1567 1568 Options Database Keys: 1569 + -options_monitor <optional filename> - prints the names and values of all runtime options as they are set. The monitor functionality is not 1570 available for options set through a file, environment variable, or on 1571 the command line. Only options set after PetscInitialize() completes will 1572 be monitored. 1573 . -options_monitor_cancel - cancel all options database monitors 1574 1575 Notes: 1576 To see all options, run your program with the -help option or consult Users-Manual: sec_gettingstarted 1577 1578 Level: intermediate 1579 1580 .keywords: set, options, database 1581 @*/ 1582 PetscErrorCode PetscOptionsSetFromOptions(PetscOptions options) 1583 { 1584 PetscBool flgc = PETSC_FALSE,flgm; 1585 PetscErrorCode ierr; 1586 char monfilename[PETSC_MAX_PATH_LEN]; 1587 PetscViewer monviewer; 1588 1589 PetscFunctionBegin; 1590 /* 1591 The options argument is currently ignored since we currently maintain only a single options database 1592 1593 options = options ? options : defaultoptions; 1594 */ 1595 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for handling options","PetscOptions");CHKERRQ(ierr); 1596 ierr = PetscOptionsString("-options_monitor","Monitor options database","PetscOptionsMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flgm);CHKERRQ(ierr); 1597 ierr = PetscOptionsBool("-options_monitor_cancel","Cancel all options database monitors","PetscOptionsMonitorCancel",flgc,&flgc,NULL);CHKERRQ(ierr); 1598 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1599 if (flgm) { 1600 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,monfilename,&monviewer);CHKERRQ(ierr); 1601 ierr = PetscOptionsMonitorSet(PetscOptionsMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 1602 } 1603 if (flgc) { ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr); } 1604 PetscFunctionReturn(0); 1605 } 1606 1607 /*@C 1608 PetscOptionsMonitorDefault - Print all options set value events. 1609 1610 Logically Collective on PETSC_COMM_WORLD 1611 1612 Input Parameters: 1613 + name - option name string 1614 . value - option value string 1615 - ctx - an ASCII viewer 1616 1617 Level: intermediate 1618 1619 .keywords: PetscOptions, default, monitor 1620 1621 .seealso: PetscOptionsMonitorSet() 1622 @*/ 1623 PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx) 1624 { 1625 PetscErrorCode ierr; 1626 PetscViewer viewer = (PetscViewer)ctx; 1627 1628 PetscFunctionBegin; 1629 if (!value) { 1630 ierr = PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);CHKERRQ(ierr); 1631 } else if (!value[0]) { 1632 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);CHKERRQ(ierr); 1633 } else { 1634 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);CHKERRQ(ierr); 1635 } 1636 PetscFunctionReturn(0); 1637 } 1638 1639 /*@C 1640 PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that 1641 modified the PETSc options database. 1642 1643 Not Collective 1644 1645 Input Parameters: 1646 + monitor - pointer to function (if this is NULL, it turns off monitoring 1647 . mctx - [optional] context for private data for the 1648 monitor routine (use NULL if no context is desired) 1649 - monitordestroy - [optional] routine that frees monitor context 1650 (may be NULL) 1651 1652 Calling Sequence of monitor: 1653 $ monitor (const char name[], const char value[], void *mctx) 1654 1655 + name - option name string 1656 . value - option value string 1657 - mctx - optional monitoring context, as set by PetscOptionsMonitorSet() 1658 1659 Options Database Keys: 1660 + -options_monitor - sets PetscOptionsMonitorDefault() 1661 - -options_monitor_cancel - cancels all monitors that have 1662 been hardwired into a code by 1663 calls to PetscOptionsMonitorSet(), but 1664 does not cancel those set via 1665 the options database. 1666 1667 Notes: 1668 The default is to do nothing. To print the name and value of options 1669 being inserted into the database, use PetscOptionsMonitorDefault() as the monitoring routine, 1670 with a null monitoring context. 1671 1672 Several different monitoring routines may be set by calling 1673 PetscOptionsMonitorSet() multiple times; all will be called in the 1674 order in which they were set. 1675 1676 Level: beginner 1677 1678 .keywords: PetscOptions, set, monitor 1679 1680 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorCancel() 1681 @*/ 1682 PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 1683 { 1684 PetscOptions options = defaultoptions; 1685 1686 PetscFunctionBegin; 1687 if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set"); 1688 options->monitor[options->numbermonitors] = monitor; 1689 options->monitordestroy[options->numbermonitors] = monitordestroy; 1690 options->monitorcontext[options->numbermonitors++] = (void*)mctx; 1691 PetscFunctionReturn(0); 1692 } 1693 1694 /*@ 1695 PetscOptionsMonitorCancel - Clears all monitors for a PetscOptions object. 1696 1697 Not Collective 1698 1699 Options Database Key: 1700 . -options_monitor_cancel - Cancels all monitors that have 1701 been hardwired into a code by calls to PetscOptionsMonitorSet(), 1702 but does not cancel those set via the options database. 1703 1704 Level: intermediate 1705 1706 .keywords: PetscOptions, set, monitor 1707 1708 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorSet() 1709 @*/ 1710 PetscErrorCode PetscOptionsMonitorCancel(void) 1711 { 1712 PetscErrorCode ierr; 1713 PetscInt i; 1714 PetscOptions options = defaultoptions; 1715 1716 PetscFunctionBegin; 1717 for (i=0; i<options->numbermonitors; i++) { 1718 if (options->monitordestroy[i]) { 1719 ierr = (*options->monitordestroy[i])(&options->monitorcontext[i]);CHKERRQ(ierr); 1720 } 1721 } 1722 options->numbermonitors = 0; 1723 PetscFunctionReturn(0); 1724 } 1725 1726 /* 1727 PetscOptionsStringToBool - Converts string to PetscBool , handles cases like "yes", "no", "true", "false", "0", "1", "off", "on". 1728 */ 1729 PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a) 1730 { 1731 PetscBool istrue,isfalse; 1732 size_t len; 1733 PetscErrorCode ierr; 1734 1735 PetscFunctionBegin; 1736 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1737 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Character string of length zero has no logical value"); 1738 ierr = PetscStrcasecmp(value,"TRUE",&istrue);CHKERRQ(ierr); 1739 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1740 ierr = PetscStrcasecmp(value,"YES",&istrue);CHKERRQ(ierr); 1741 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1742 ierr = PetscStrcasecmp(value,"1",&istrue);CHKERRQ(ierr); 1743 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1744 ierr = PetscStrcasecmp(value,"on",&istrue);CHKERRQ(ierr); 1745 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1746 ierr = PetscStrcasecmp(value,"FALSE",&isfalse);CHKERRQ(ierr); 1747 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1748 ierr = PetscStrcasecmp(value,"NO",&isfalse);CHKERRQ(ierr); 1749 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1750 ierr = PetscStrcasecmp(value,"0",&isfalse);CHKERRQ(ierr); 1751 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1752 ierr = PetscStrcasecmp(value,"off",&isfalse);CHKERRQ(ierr); 1753 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1754 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value); 1755 } 1756 1757 /* 1758 PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide" 1759 */ 1760 PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a) 1761 { 1762 PetscErrorCode ierr; 1763 size_t len; 1764 PetscBool decide,tdefault,mouse; 1765 1766 PetscFunctionBegin; 1767 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1768 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 1769 1770 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);CHKERRQ(ierr); 1771 if (!tdefault) { 1772 ierr = PetscStrcasecmp(name,"DEFAULT",&tdefault);CHKERRQ(ierr); 1773 } 1774 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&decide);CHKERRQ(ierr); 1775 if (!decide) { 1776 ierr = PetscStrcasecmp(name,"DECIDE",&decide);CHKERRQ(ierr); 1777 } 1778 ierr = PetscStrcasecmp(name,"mouse",&mouse);CHKERRQ(ierr); 1779 1780 if (tdefault) *a = PETSC_DEFAULT; 1781 else if (decide) *a = PETSC_DECIDE; 1782 else if (mouse) *a = -1; 1783 else { 1784 char *endptr; 1785 long strtolval; 1786 1787 strtolval = strtol(name,&endptr,10); 1788 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); 1789 1790 #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL) 1791 (void) strtolval; 1792 *a = atoll(name); 1793 #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64) 1794 (void) strtolval; 1795 *a = _atoi64(name); 1796 #else 1797 *a = (PetscInt)strtolval; 1798 #endif 1799 } 1800 PetscFunctionReturn(0); 1801 } 1802 1803 #if defined(PETSC_USE_REAL___FLOAT128) 1804 #include <quadmath.h> 1805 #endif 1806 1807 static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr) 1808 { 1809 PetscFunctionBegin; 1810 #if defined(PETSC_USE_REAL___FLOAT128) 1811 *a = strtoflt128(name,endptr); 1812 #else 1813 *a = (PetscReal)strtod(name,endptr); 1814 #endif 1815 PetscFunctionReturn(0); 1816 } 1817 1818 static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary) 1819 { 1820 PetscBool hasi = PETSC_FALSE; 1821 char *ptr; 1822 PetscReal strtoval; 1823 PetscErrorCode ierr; 1824 1825 PetscFunctionBegin; 1826 ierr = PetscStrtod(name,&strtoval,&ptr);CHKERRQ(ierr); 1827 if (ptr == name) { 1828 strtoval = 1.; 1829 hasi = PETSC_TRUE; 1830 if (name[0] == 'i') { 1831 ptr++; 1832 } else if (name[0] == '+' && name[1] == 'i') { 1833 ptr += 2; 1834 } else if (name[0] == '-' && name[1] == 'i') { 1835 strtoval = -1.; 1836 ptr += 2; 1837 } 1838 } else if (*ptr == 'i') { 1839 hasi = PETSC_TRUE; 1840 ptr++; 1841 } 1842 *endptr = ptr; 1843 *isImaginary = hasi; 1844 if (hasi) { 1845 #if !defined(PETSC_USE_COMPLEX) 1846 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name); 1847 #else 1848 *a = PetscCMPLX(0.,strtoval); 1849 #endif 1850 } else { 1851 *a = strtoval; 1852 } 1853 PetscFunctionReturn(0); 1854 } 1855 1856 /* 1857 Converts a string to PetscReal value. Handles special cases like "default" and "decide" 1858 */ 1859 PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a) 1860 { 1861 size_t len; 1862 PetscBool match; 1863 char *endptr; 1864 PetscErrorCode ierr; 1865 1866 PetscFunctionBegin; 1867 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1868 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"String of length zero has no numerical value"); 1869 1870 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&match);CHKERRQ(ierr); 1871 if (!match) { 1872 ierr = PetscStrcasecmp(name,"DEFAULT",&match);CHKERRQ(ierr); 1873 } 1874 if (match) {*a = PETSC_DEFAULT; PetscFunctionReturn(0);} 1875 1876 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&match);CHKERRQ(ierr); 1877 if (!match) { 1878 ierr = PetscStrcasecmp(name,"DECIDE",&match);CHKERRQ(ierr); 1879 } 1880 if (match) {*a = PETSC_DECIDE; PetscFunctionReturn(0);} 1881 1882 ierr = PetscStrtod(name,a,&endptr);CHKERRQ(ierr); 1883 if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name); 1884 PetscFunctionReturn(0); 1885 } 1886 1887 PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a) 1888 { 1889 PetscBool imag1; 1890 size_t len; 1891 PetscScalar val = 0.; 1892 char *ptr = NULL; 1893 PetscErrorCode ierr; 1894 1895 PetscFunctionBegin; 1896 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1897 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 1898 ierr = PetscStrtoz(name,&val,&ptr,&imag1);CHKERRQ(ierr); 1899 #if defined(PETSC_USE_COMPLEX) 1900 if ((size_t) (ptr - name) < len) { 1901 PetscBool imag2; 1902 PetscScalar val2; 1903 1904 ierr = PetscStrtoz(ptr,&val2,&ptr,&imag2);CHKERRQ(ierr); 1905 if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name); 1906 val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2)); 1907 } 1908 #endif 1909 if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name); 1910 *a = val; 1911 PetscFunctionReturn(0); 1912 } 1913 1914 /*@C 1915 PetscOptionsGetBool - Gets the Logical (true or false) value for a particular 1916 option in the database. 1917 1918 Not Collective 1919 1920 Input Parameters: 1921 + options - options database, use NULL for default global database 1922 . pre - the string to prepend to the name or NULL 1923 - name - the option one is seeking 1924 1925 Output Parameter: 1926 + ivalue - the logical value to return 1927 - set - PETSC_TRUE if found, else PETSC_FALSE 1928 1929 Level: beginner 1930 1931 Notes: 1932 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 1933 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 1934 1935 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 1936 is equivalent to -requested_bool true 1937 1938 If the user does not supply the option at all ivalue is NOT changed. Thus 1939 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 1940 1941 Concepts: options database^has logical 1942 1943 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1944 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(), 1945 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1946 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1947 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1948 PetscOptionsFList(), PetscOptionsEList() 1949 @*/ 1950 PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set) 1951 { 1952 const char *value; 1953 PetscBool flag; 1954 PetscErrorCode ierr; 1955 1956 PetscFunctionBegin; 1957 PetscValidCharPointer(name,3); 1958 if (ivalue) PetscValidIntPointer(ivalue,4); 1959 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 1960 if (flag) { 1961 if (set) *set = PETSC_TRUE; 1962 if (!value) { 1963 if (ivalue) *ivalue = PETSC_TRUE; 1964 } else { 1965 ierr = PetscOptionsStringToBool(value, &flag);CHKERRQ(ierr); 1966 if (ivalue) *ivalue = flag; 1967 } 1968 } else { 1969 if (set) *set = PETSC_FALSE; 1970 } 1971 PetscFunctionReturn(0); 1972 } 1973 1974 /*@C 1975 PetscOptionsGetEList - Puts a list of option values that a single one may be selected from 1976 1977 Not Collective 1978 1979 Input Parameters: 1980 + options - options database, use NULL for default global database 1981 . pre - the string to prepend to the name or NULL 1982 . opt - option name 1983 . list - the possible choices (one of these must be selected, anything else is invalid) 1984 . ntext - number of choices 1985 1986 Output Parameter: 1987 + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed) 1988 - set - PETSC_TRUE if found, else PETSC_FALSE 1989 1990 Level: intermediate 1991 1992 Notes: 1993 If the user does not supply the option value is NOT changed. Thus 1994 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 1995 1996 See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 1997 1998 Concepts: options database^list 1999 2000 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2001 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2002 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2003 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2004 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2005 PetscOptionsFList(), PetscOptionsEList() 2006 @*/ 2007 PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set) 2008 { 2009 PetscErrorCode ierr; 2010 size_t alen,len = 0; 2011 char *svalue; 2012 PetscBool aset,flg = PETSC_FALSE; 2013 PetscInt i; 2014 2015 PetscFunctionBegin; 2016 PetscValidCharPointer(opt,3); 2017 for (i=0; i<ntext; i++) { 2018 ierr = PetscStrlen(list[i],&alen);CHKERRQ(ierr); 2019 if (alen > len) len = alen; 2020 } 2021 len += 5; /* a little extra space for user mistypes */ 2022 ierr = PetscMalloc1(len,&svalue);CHKERRQ(ierr); 2023 ierr = PetscOptionsGetString(options,pre,opt,svalue,len,&aset);CHKERRQ(ierr); 2024 if (aset) { 2025 ierr = PetscEListFind(ntext,list,svalue,value,&flg);CHKERRQ(ierr); 2026 if (!flg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s",svalue,pre ? pre : "",opt+1); 2027 if (set) *set = PETSC_TRUE; 2028 } else if (set) *set = PETSC_FALSE; 2029 ierr = PetscFree(svalue);CHKERRQ(ierr); 2030 PetscFunctionReturn(0); 2031 } 2032 2033 /*@C 2034 PetscOptionsGetEnum - Gets the enum value for a particular option in the database. 2035 2036 Not Collective 2037 2038 Input Parameters: 2039 + options - options database, use NULL for default global database 2040 . pre - option prefix or NULL 2041 . opt - option name 2042 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2043 - defaultv - the default (current) value 2044 2045 Output Parameter: 2046 + value - the value to return 2047 - set - PETSC_TRUE if found, else PETSC_FALSE 2048 2049 Level: beginner 2050 2051 Concepts: options database 2052 2053 Notes: 2054 If the user does not supply the option value is NOT changed. Thus 2055 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2056 2057 List is usually something like PCASMTypes or some other predefined list of enum names 2058 2059 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2060 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2061 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2062 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2063 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2064 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2065 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2066 @*/ 2067 PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set) 2068 { 2069 PetscErrorCode ierr; 2070 PetscInt ntext = 0,tval; 2071 PetscBool fset; 2072 2073 PetscFunctionBegin; 2074 PetscValidCharPointer(opt,3); 2075 while (list[ntext++]) { 2076 if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 2077 } 2078 if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 2079 ntext -= 3; 2080 ierr = PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);CHKERRQ(ierr); 2081 /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 2082 if (fset) *value = (PetscEnum)tval; 2083 if (set) *set = fset; 2084 PetscFunctionReturn(0); 2085 } 2086 2087 /*@C 2088 PetscOptionsGetInt - Gets the integer value for a particular option in the database. 2089 2090 Not Collective 2091 2092 Input Parameters: 2093 + options - options database, use NULL for default global database 2094 . pre - the string to prepend to the name or NULL 2095 - name - the option one is seeking 2096 2097 Output Parameter: 2098 + ivalue - the integer value to return 2099 - set - PETSC_TRUE if found, else PETSC_FALSE 2100 2101 Level: beginner 2102 2103 Notes: 2104 If the user does not supply the option ivalue is NOT changed. Thus 2105 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2106 2107 Concepts: options database^has int 2108 2109 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 2110 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2111 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2112 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2113 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2114 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2115 PetscOptionsFList(), PetscOptionsEList() 2116 @*/ 2117 PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set) 2118 { 2119 const char *value; 2120 PetscErrorCode ierr; 2121 PetscBool flag; 2122 2123 PetscFunctionBegin; 2124 PetscValidCharPointer(name,3); 2125 PetscValidIntPointer(ivalue,4); 2126 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2127 if (flag) { 2128 if (!value) { 2129 if (set) *set = PETSC_FALSE; 2130 } else { 2131 if (set) *set = PETSC_TRUE; 2132 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2133 } 2134 } else { 2135 if (set) *set = PETSC_FALSE; 2136 } 2137 PetscFunctionReturn(0); 2138 } 2139 2140 /*@C 2141 PetscOptionsGetReal - Gets the double precision value for a particular 2142 option in the database. 2143 2144 Not Collective 2145 2146 Input Parameters: 2147 + options - options database, use NULL for default global database 2148 . pre - string to prepend to each name or NULL 2149 - name - the option one is seeking 2150 2151 Output Parameter: 2152 + dvalue - the double value to return 2153 - set - PETSC_TRUE if found, PETSC_FALSE if not found 2154 2155 Notes: 2156 If the user does not supply the option dvalue is NOT changed. Thus 2157 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2158 2159 Level: beginner 2160 2161 Concepts: options database^has double 2162 2163 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2164 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(), 2165 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2166 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2167 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2168 PetscOptionsFList(), PetscOptionsEList() 2169 @*/ 2170 PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set) 2171 { 2172 const char *value; 2173 PetscBool flag; 2174 PetscErrorCode ierr; 2175 2176 PetscFunctionBegin; 2177 PetscValidCharPointer(name,3); 2178 PetscValidRealPointer(dvalue,4); 2179 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2180 if (flag) { 2181 if (!value) { 2182 if (set) *set = PETSC_FALSE; 2183 } else { 2184 if (set) *set = PETSC_TRUE; 2185 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2186 } 2187 } else { 2188 if (set) *set = PETSC_FALSE; 2189 } 2190 PetscFunctionReturn(0); 2191 } 2192 2193 /*@C 2194 PetscOptionsGetScalar - Gets the scalar value for a particular 2195 option in the database. 2196 2197 Not Collective 2198 2199 Input Parameters: 2200 + options - options database, use NULL for default global database 2201 . pre - string to prepend to each name or NULL 2202 - name - the option one is seeking 2203 2204 Output Parameter: 2205 + dvalue - the double value to return 2206 - set - PETSC_TRUE if found, else PETSC_FALSE 2207 2208 Level: beginner 2209 2210 Usage: 2211 A complex number 2+3i must be specified with NO spaces 2212 2213 Notes: 2214 If the user does not supply the option dvalue is NOT changed. Thus 2215 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2216 2217 Concepts: options database^has scalar 2218 2219 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2220 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2221 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2222 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2223 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2224 PetscOptionsFList(), PetscOptionsEList() 2225 @*/ 2226 PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set) 2227 { 2228 const char *value; 2229 PetscBool flag; 2230 PetscErrorCode ierr; 2231 2232 PetscFunctionBegin; 2233 PetscValidCharPointer(name,3); 2234 PetscValidScalarPointer(dvalue,4); 2235 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2236 if (flag) { 2237 if (!value) { 2238 if (set) *set = PETSC_FALSE; 2239 } else { 2240 #if !defined(PETSC_USE_COMPLEX) 2241 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2242 #else 2243 ierr = PetscOptionsStringToScalar(value,dvalue);CHKERRQ(ierr); 2244 #endif 2245 if (set) *set = PETSC_TRUE; 2246 } 2247 } else { /* flag */ 2248 if (set) *set = PETSC_FALSE; 2249 } 2250 PetscFunctionReturn(0); 2251 } 2252 2253 /*@C 2254 PetscOptionsGetString - Gets the string value for a particular option in 2255 the database. 2256 2257 Not Collective 2258 2259 Input Parameters: 2260 + options - options database, use NULL for default global database 2261 . pre - string to prepend to name or NULL 2262 . name - the option one is seeking 2263 - len - maximum length of the string including null termination 2264 2265 Output Parameters: 2266 + string - location to copy string 2267 - set - PETSC_TRUE if found, else PETSC_FALSE 2268 2269 Level: beginner 2270 2271 Fortran Note: 2272 The Fortran interface is slightly different from the C/C++ 2273 interface (len is not used). Sample usage in Fortran follows 2274 .vb 2275 character *20 string 2276 PetscErrorCode ierr 2277 PetscBool set 2278 call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr) 2279 .ve 2280 2281 Notes: 2282 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 2283 2284 If the user does not use the option then the string is not changed. Thus 2285 you should ALWAYS initialize the string if you access it without first checking if the set flag is true. 2286 2287 Concepts: options database^string 2288 2289 Note: 2290 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). 2291 2292 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2293 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2294 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2295 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2296 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2297 PetscOptionsFList(), PetscOptionsEList() 2298 @*/ 2299 PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set) 2300 { 2301 const char *value; 2302 PetscBool flag; 2303 PetscErrorCode ierr; 2304 2305 PetscFunctionBegin; 2306 PetscValidCharPointer(name,3); 2307 PetscValidCharPointer(string,4); 2308 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2309 if (!flag) { 2310 if (set) *set = PETSC_FALSE; 2311 } else { 2312 if (set) *set = PETSC_TRUE; 2313 if (value) { 2314 ierr = PetscStrncpy(string,value,len);CHKERRQ(ierr); 2315 } else { 2316 ierr = PetscMemzero(string,len);CHKERRQ(ierr); 2317 } 2318 } 2319 PetscFunctionReturn(0); 2320 } 2321 2322 char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[]) 2323 { 2324 const char *value; 2325 PetscBool flag; 2326 PetscErrorCode ierr; 2327 2328 PetscFunctionBegin; 2329 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(0); 2330 if (flag) PetscFunctionReturn((char*)value); 2331 else PetscFunctionReturn(0); 2332 } 2333 2334 /*@C 2335 PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular 2336 option in the database. The values must be separated with commas with 2337 no intervening spaces. 2338 2339 Not Collective 2340 2341 Input Parameters: 2342 + options - options database, use NULL for default global database 2343 . pre - string to prepend to each name or NULL 2344 . name - the option one is seeking 2345 - nmax - maximum number of values to retrieve 2346 2347 Output Parameter: 2348 + dvalue - the integer values to return 2349 . nmax - actual number of values retreived 2350 - set - PETSC_TRUE if found, else PETSC_FALSE 2351 2352 Level: beginner 2353 2354 Concepts: options database^array of ints 2355 2356 Notes: 2357 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 2358 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 2359 2360 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2361 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2362 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2363 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2364 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2365 PetscOptionsFList(), PetscOptionsEList() 2366 @*/ 2367 PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set) 2368 { 2369 const char *svalue; 2370 char *value; 2371 PetscErrorCode ierr; 2372 PetscInt n = 0; 2373 PetscBool flag; 2374 PetscToken token; 2375 2376 PetscFunctionBegin; 2377 PetscValidCharPointer(name,3); 2378 PetscValidIntPointer(dvalue,4); 2379 PetscValidIntPointer(nmax,5); 2380 2381 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2382 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2383 if (set) *set = PETSC_TRUE; 2384 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2385 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2386 while (value && n < *nmax) { 2387 ierr = PetscOptionsStringToBool(value,dvalue);CHKERRQ(ierr); 2388 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2389 dvalue++; 2390 n++; 2391 } 2392 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2393 *nmax = n; 2394 PetscFunctionReturn(0); 2395 } 2396 2397 /*@C 2398 PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database. 2399 2400 Not Collective 2401 2402 Input Parameters: 2403 + options - options database, use NULL for default global database 2404 . pre - option prefix or NULL 2405 . name - option name 2406 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2407 - nmax - maximum number of values to retrieve 2408 2409 Output Parameters: 2410 + ivalue - the enum values to return 2411 . nmax - actual number of values retreived 2412 - set - PETSC_TRUE if found, else PETSC_FALSE 2413 2414 Level: beginner 2415 2416 Concepts: options database 2417 2418 Notes: 2419 The array must be passed as a comma separated list. 2420 2421 There must be no intervening spaces between the values. 2422 2423 list is usually something like PCASMTypes or some other predefined list of enum names. 2424 2425 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2426 PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2427 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(), 2428 PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(), 2429 PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2430 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2431 @*/ 2432 PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set) 2433 { 2434 const char *svalue; 2435 char *value; 2436 PetscInt n = 0; 2437 PetscEnum evalue; 2438 PetscBool flag; 2439 PetscToken token; 2440 PetscErrorCode ierr; 2441 2442 PetscFunctionBegin; 2443 PetscValidCharPointer(name,3); 2444 PetscValidPointer(list,4); 2445 PetscValidPointer(ivalue,5); 2446 PetscValidIntPointer(nmax,6); 2447 2448 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2449 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2450 if (set) *set = PETSC_TRUE; 2451 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2452 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2453 while (value && n < *nmax) { 2454 ierr = PetscEnumFind(list,value,&evalue,&flag);CHKERRQ(ierr); 2455 if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1); 2456 ivalue[n++] = evalue; 2457 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2458 } 2459 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2460 *nmax = n; 2461 PetscFunctionReturn(0); 2462 } 2463 2464 /*@C 2465 PetscOptionsGetIntArray - Gets an array of integer values for a particular 2466 option in the database. 2467 2468 Not Collective 2469 2470 Input Parameters: 2471 + options - options database, use NULL for default global database 2472 . pre - string to prepend to each name or NULL 2473 . name - the option one is seeking 2474 - nmax - maximum number of values to retrieve 2475 2476 Output Parameter: 2477 + ivalue - the integer values to return 2478 . nmax - actual number of values retreived 2479 - set - PETSC_TRUE if found, else PETSC_FALSE 2480 2481 Level: beginner 2482 2483 Notes: 2484 The array can be passed as 2485 a comma separated list: 0,1,2,3,4,5,6,7 2486 a range (start-end+1): 0-8 2487 a range with given increment (start-end+1:inc): 0-7:2 2488 a combination of values and ranges separated by commas: 0,1-8,8-15:2 2489 2490 There must be no intervening spaces between the values. 2491 2492 Concepts: options database^array of ints 2493 2494 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2495 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2496 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2497 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2498 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2499 PetscOptionsFList(), PetscOptionsEList() 2500 @*/ 2501 PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set) 2502 { 2503 const char *svalue; 2504 char *value; 2505 PetscErrorCode ierr; 2506 PetscInt n = 0,i,j,start,end,inc,nvalues; 2507 size_t len; 2508 PetscBool flag,foundrange; 2509 PetscToken token; 2510 2511 PetscFunctionBegin; 2512 PetscValidCharPointer(name,3); 2513 PetscValidIntPointer(ivalue,4); 2514 PetscValidIntPointer(nmax,5); 2515 2516 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2517 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2518 if (set) *set = PETSC_TRUE; 2519 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2520 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2521 while (value && n < *nmax) { 2522 /* look for form d-D where d and D are integers */ 2523 foundrange = PETSC_FALSE; 2524 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 2525 if (value[0] == '-') i=2; 2526 else i=1; 2527 for (;i<(int)len; i++) { 2528 if (value[i] == '-') { 2529 if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 2530 value[i] = 0; 2531 2532 ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 2533 inc = 1; 2534 j = i+1; 2535 for (;j<(int)len; j++) { 2536 if (value[j] == ':') { 2537 value[j] = 0; 2538 2539 ierr = PetscOptionsStringToInt(value+j+1,&inc);CHKERRQ(ierr); 2540 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); 2541 break; 2542 } 2543 } 2544 ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 2545 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); 2546 nvalues = (end-start)/inc + (end-start)%inc; 2547 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); 2548 for (;start<end; start+=inc) { 2549 *ivalue = start; ivalue++;n++; 2550 } 2551 foundrange = PETSC_TRUE; 2552 break; 2553 } 2554 } 2555 if (!foundrange) { 2556 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2557 ivalue++; 2558 n++; 2559 } 2560 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2561 } 2562 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2563 *nmax = n; 2564 PetscFunctionReturn(0); 2565 } 2566 2567 /*@C 2568 PetscOptionsGetRealArray - Gets an array of double precision values for a 2569 particular option in the database. The values must be separated with 2570 commas with no intervening spaces. 2571 2572 Not Collective 2573 2574 Input Parameters: 2575 + options - options database, use NULL for default global database 2576 . pre - string to prepend to each name or NULL 2577 . name - the option one is seeking 2578 - nmax - maximum number of values to retrieve 2579 2580 Output Parameters: 2581 + dvalue - the double values to return 2582 . nmax - actual number of values retreived 2583 - set - PETSC_TRUE if found, else PETSC_FALSE 2584 2585 Level: beginner 2586 2587 Concepts: options database^array of doubles 2588 2589 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2590 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2591 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2592 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2593 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2594 PetscOptionsFList(), PetscOptionsEList() 2595 @*/ 2596 PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set) 2597 { 2598 const char *svalue; 2599 char *value; 2600 PetscErrorCode ierr; 2601 PetscInt n = 0; 2602 PetscBool flag; 2603 PetscToken token; 2604 2605 PetscFunctionBegin; 2606 PetscValidCharPointer(name,3); 2607 PetscValidRealPointer(dvalue,4); 2608 PetscValidIntPointer(nmax,5); 2609 2610 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2611 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2612 if (set) *set = PETSC_TRUE; 2613 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2614 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2615 while (value && n < *nmax) { 2616 ierr = PetscOptionsStringToReal(value,dvalue++);CHKERRQ(ierr); 2617 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2618 n++; 2619 } 2620 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2621 *nmax = n; 2622 PetscFunctionReturn(0); 2623 } 2624 2625 /*@C 2626 PetscOptionsGetScalarArray - Gets an array of scalars for a 2627 particular option in the database. The values must be separated with 2628 commas with no intervening spaces. 2629 2630 Not Collective 2631 2632 Input Parameters: 2633 + options - options database, use NULL for default global database 2634 . pre - string to prepend to each name or NULL 2635 . name - the option one is seeking 2636 - nmax - maximum number of values to retrieve 2637 2638 Output Parameters: 2639 + dvalue - the scalar values to return 2640 . nmax - actual number of values retreived 2641 - set - PETSC_TRUE if found, else PETSC_FALSE 2642 2643 Level: beginner 2644 2645 Concepts: options database^array of doubles 2646 2647 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2648 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2649 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2650 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2651 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2652 PetscOptionsFList(), PetscOptionsEList() 2653 @*/ 2654 PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],PetscInt *nmax,PetscBool *set) 2655 { 2656 const char *svalue; 2657 char *value; 2658 PetscErrorCode ierr; 2659 PetscInt n = 0; 2660 PetscBool flag; 2661 PetscToken token; 2662 2663 PetscFunctionBegin; 2664 PetscValidCharPointer(name,3); 2665 PetscValidRealPointer(dvalue,4); 2666 PetscValidIntPointer(nmax,5); 2667 2668 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2669 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2670 if (set) *set = PETSC_TRUE; 2671 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2672 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2673 while (value && n < *nmax) { 2674 ierr = PetscOptionsStringToScalar(value,dvalue++);CHKERRQ(ierr); 2675 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2676 n++; 2677 } 2678 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2679 *nmax = n; 2680 PetscFunctionReturn(0); 2681 } 2682 2683 /*@C 2684 PetscOptionsGetStringArray - Gets an array of string values for a particular 2685 option in the database. The values must be separated with commas with 2686 no intervening spaces. 2687 2688 Not Collective 2689 2690 Input Parameters: 2691 + options - options database, use NULL for default global database 2692 . pre - string to prepend to name or NULL 2693 . name - the option one is seeking 2694 - nmax - maximum number of strings 2695 2696 Output Parameter: 2697 + strings - location to copy strings 2698 - set - PETSC_TRUE if found, else PETSC_FALSE 2699 2700 Level: beginner 2701 2702 Notes: 2703 The user should pass in an array of pointers to char, to hold all the 2704 strings returned by this function. 2705 2706 The user is responsible for deallocating the strings that are 2707 returned. The Fortran interface for this routine is not supported. 2708 2709 Contributed by Matthew Knepley. 2710 2711 Concepts: options database^array of strings 2712 2713 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2714 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2715 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2716 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2717 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2718 PetscOptionsFList(), PetscOptionsEList() 2719 @*/ 2720 PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set) 2721 { 2722 const char *svalue; 2723 char *value; 2724 PetscErrorCode ierr; 2725 PetscInt n = 0; 2726 PetscBool flag; 2727 PetscToken token; 2728 2729 PetscFunctionBegin; 2730 PetscValidCharPointer(name,3); 2731 PetscValidPointer(strings,4); 2732 PetscValidIntPointer(nmax,5); 2733 2734 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2735 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2736 if (set) *set = PETSC_TRUE; 2737 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2738 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2739 while (value && n < *nmax) { 2740 ierr = PetscStrallocpy(value,&strings[n]);CHKERRQ(ierr); 2741 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2742 n++; 2743 } 2744 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2745 *nmax = n; 2746 PetscFunctionReturn(0); 2747 } 2748 2749 /*@C 2750 PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one 2751 2752 Prints a deprecation warning, unless an option is supplied to suppress. 2753 2754 Not Collective 2755 2756 Input Parameters: 2757 + pre - string to prepend to name or NULL 2758 . oldname - the old, deprecated option 2759 . newname - the new option, or NULL if option is purely removed 2760 . version - a string describing the version of first deprecation, e.g. "3.9" 2761 - info - additional information string, or NULL. 2762 2763 Options Database Keys: 2764 . -options_suppress_deprecated_warnings - do not print deprecation warnings 2765 2766 Notes: 2767 Must be called between PetscOptionsBegin() and PetscOptionsEnd(). 2768 If newname is provided, the old option is replaced. Otherwise, it remains 2769 in the options database. 2770 If an option is not replaced, the info argument should be used to advise the user 2771 on how to proceed. 2772 There is a limit on the length of the warning printed, so very long strings 2773 provided as info may be truncated. 2774 2775 Level: developer 2776 2777 .seealso: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsScalar(), PetscOptionsBool(), PetscOptionsString(), PetscOptionsSetValue() 2778 2779 @*/ 2780 PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[]) 2781 { 2782 PetscErrorCode ierr; 2783 PetscBool found,quiet; 2784 const char *value; 2785 const char * const quietopt="-options_suppress_deprecated_warnings"; 2786 char msg[4096]; 2787 2788 PetscFunctionBegin; 2789 PetscValidCharPointer(oldname,2); 2790 PetscValidCharPointer(version,4); 2791 2792 ierr = PetscOptionsFindPair(PetscOptionsObject->options,PetscOptionsObject->prefix,oldname,&value,&found);CHKERRQ(ierr); 2793 if (found) { 2794 if (newname) { 2795 if (PetscOptionsObject->prefix) { 2796 ierr = PetscOptionsPrefixPush(PetscOptionsObject->options,PetscOptionsObject->prefix);CHKERRQ(ierr); 2797 } 2798 ierr = PetscOptionsSetValue(PetscOptionsObject->options,newname,value);CHKERRQ(ierr); 2799 if (PetscOptionsObject->prefix) { 2800 ierr = PetscOptionsPrefixPop(PetscOptionsObject->options);CHKERRQ(ierr); 2801 } 2802 ierr = PetscOptionsClearValue(PetscOptionsObject->options,oldname);CHKERRQ(ierr); 2803 } 2804 quiet = PETSC_FALSE; 2805 ierr = PetscOptionsGetBool(PetscOptionsObject->options,NULL,quietopt,&quiet,NULL);CHKERRQ(ierr); 2806 if (!quiet) { 2807 ierr = PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");CHKERRQ(ierr); 2808 ierr = PetscStrcat(msg,oldname);CHKERRQ(ierr); 2809 ierr = PetscStrcat(msg," is deprecated as of version ");CHKERRQ(ierr); 2810 ierr = PetscStrcat(msg,version);CHKERRQ(ierr); 2811 ierr = PetscStrcat(msg," and will be removed in a future release.");CHKERRQ(ierr); 2812 if (newname) { 2813 ierr = PetscStrcat(msg," Please use the option ");CHKERRQ(ierr); 2814 ierr = PetscStrcat(msg,newname);CHKERRQ(ierr); 2815 ierr = PetscStrcat(msg," instead.");CHKERRQ(ierr); 2816 } 2817 if (info) { 2818 ierr = PetscStrcat(msg," ");CHKERRQ(ierr); 2819 ierr = PetscStrcat(msg,info);CHKERRQ(ierr); 2820 } 2821 ierr = PetscStrcat(msg," (Silence this warning with ");CHKERRQ(ierr); 2822 ierr = PetscStrcat(msg,quietopt);CHKERRQ(ierr); 2823 ierr = PetscStrcat(msg,")\n");CHKERRQ(ierr); 2824 ierr = PetscPrintf(PetscOptionsObject->comm,msg);CHKERRQ(ierr); 2825 } 2826 } 2827 PetscFunctionReturn(0); 2828 } 2829