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