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