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 *eoptions = NULL; 701 size_t len = 0; 702 if (!rank) { 703 eoptions = (char*)getenv("PETSC_OPTIONS_YAML"); 704 ierr = PetscStrlen(eoptions,&len);CHKERRQ(ierr); 705 ierr = MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 706 } else { 707 ierr = MPI_Bcast(&len,1,MPIU_SIZE_T,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 708 if (len) { 709 ierr = PetscMalloc1(len+1,&eoptions);CHKERRQ(ierr); 710 } 711 } 712 if (len) { 713 ierr = MPI_Bcast(eoptions,len,MPI_CHAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 714 if (rank) eoptions[len] = 0; 715 ierr = PetscOptionsInsertStringYAML(options,eoptions);CHKERRQ(ierr); 716 if (rank) {ierr = PetscFree(eoptions);CHKERRQ(ierr);} 717 } 718 } 719 { 720 char yaml_file[PETSC_MAX_PATH_LEN]; 721 char yaml_string[BUFSIZ]; 722 PetscBool yaml_flg; 723 ierr = PetscOptionsGetString(NULL,NULL,"-options_file_yaml",yaml_file,PETSC_MAX_PATH_LEN,&yaml_flg);CHKERRQ(ierr); 724 if (yaml_flg) { 725 ierr = PetscOptionsInsertFileYAML(PETSC_COMM_WORLD,yaml_file,PETSC_TRUE);CHKERRQ(ierr); 726 } 727 ierr = PetscOptionsGetString(NULL,NULL,"-options_string_yaml",yaml_string,BUFSIZ,&yaml_flg);CHKERRQ(ierr); 728 if (yaml_flg) { 729 ierr = PetscOptionsInsertStringYAML(NULL,yaml_string);CHKERRQ(ierr); 730 } 731 } 732 #endif 733 734 /* insert command line options again because they take precedence over arguments in petscrc/environment */ 735 if (argc && args && *argc) {ierr = PetscOptionsInsertArgs(options,*argc,*args);CHKERRQ(ierr);} 736 PetscFunctionReturn(0); 737 } 738 739 /*@C 740 PetscOptionsView - Prints the options that have been loaded. This is 741 useful for debugging purposes. 742 743 Logically Collective on PetscViewer 744 745 Input Parameter: 746 + options - options database, use NULL for default global database 747 - viewer - must be an PETSCVIEWERASCII viewer 748 749 Options Database Key: 750 . -options_view - Activates PetscOptionsView() within PetscFinalize() 751 752 Notes: 753 Only the rank zero process of MPI_Comm used to create view prints the option values. Other processes 754 may have different values but they are not printed. 755 756 Level: advanced 757 758 .seealso: PetscOptionsAllUsed() 759 @*/ 760 PetscErrorCode PetscOptionsView(PetscOptions options,PetscViewer viewer) 761 { 762 PetscErrorCode ierr; 763 PetscInt i; 764 PetscBool isascii; 765 766 PetscFunctionBegin; 767 if (viewer) PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 768 options = options ? options : defaultoptions; 769 if (!viewer) viewer = PETSC_VIEWER_STDOUT_WORLD; 770 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 771 if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Only supports ASCII viewer"); 772 773 if (!options->N) { 774 ierr = PetscViewerASCIIPrintf(viewer,"#No PETSc Option Table entries\n");CHKERRQ(ierr); 775 PetscFunctionReturn(0); 776 } 777 778 ierr = PetscViewerASCIIPrintf(viewer,"#PETSc Option Table entries:\n");CHKERRQ(ierr); 779 for (i=0; i<options->N; i++) { 780 if (options->values[i]) { 781 ierr = PetscViewerASCIIPrintf(viewer,"-%s %s\n",options->names[i],options->values[i]);CHKERRQ(ierr); 782 } else { 783 ierr = PetscViewerASCIIPrintf(viewer,"-%s\n",options->names[i]);CHKERRQ(ierr); 784 } 785 } 786 ierr = PetscViewerASCIIPrintf(viewer,"#End of PETSc Option Table entries\n");CHKERRQ(ierr); 787 PetscFunctionReturn(0); 788 } 789 790 /* 791 Called by error handlers to print options used in run 792 */ 793 PETSC_EXTERN PetscErrorCode PetscOptionsViewError(void) 794 { 795 PetscInt i; 796 PetscOptions options = defaultoptions; 797 798 PetscFunctionBegin; 799 if (options->N) { 800 (*PetscErrorPrintf)("PETSc Option Table entries:\n"); 801 } else { 802 (*PetscErrorPrintf)("No PETSc Option Table entries\n"); 803 } 804 for (i=0; i<options->N; i++) { 805 if (options->values[i]) { 806 (*PetscErrorPrintf)("-%s %s\n",options->names[i],options->values[i]); 807 } else { 808 (*PetscErrorPrintf)("-%s\n",options->names[i]); 809 } 810 } 811 PetscFunctionReturn(0); 812 } 813 814 /*@C 815 PetscOptionsPrefixPush - Designate a prefix to be used by all options insertions to follow. 816 817 Logically Collective 818 819 Input Parameter: 820 + options - options database, or NULL for the default global database 821 - prefix - The string to append to the existing prefix 822 823 Options Database Keys: 824 + -prefix_push <some_prefix_> - push the given prefix 825 - -prefix_pop - pop the last prefix 826 827 Notes: 828 It is common to use this in conjunction with -options_file as in 829 830 $ -prefix_push system1_ -options_file system1rc -prefix_pop -prefix_push system2_ -options_file system2rc -prefix_pop 831 832 where the files no longer require all options to be prefixed with -system2_. 833 834 The collectivity of this routine is complex; only the MPI processes that call this routine will 835 have the affect of these options. If some processes that create objects call this routine and others do 836 not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options 837 on different ranks. 838 839 Level: advanced 840 841 .seealso: PetscOptionsPrefixPop(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue() 842 @*/ 843 PetscErrorCode PetscOptionsPrefixPush(PetscOptions options,const char prefix[]) 844 { 845 PetscErrorCode ierr; 846 size_t n; 847 PetscInt start; 848 char key[MAXOPTNAME+1]; 849 PetscBool valid; 850 851 PetscFunctionBegin; 852 PetscValidCharPointer(prefix,1); 853 options = options ? options : defaultoptions; 854 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); 855 key[0] = '-'; /* keys must start with '-' */ 856 ierr = PetscStrncpy(key+1,prefix,sizeof(key)-1);CHKERRQ(ierr); 857 ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr); 858 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); 859 start = options->prefixind ? options->prefixstack[options->prefixind-1] : 0; 860 ierr = PetscStrlen(prefix,&n);CHKERRQ(ierr); 861 if (n+1 > sizeof(options->prefix)-start) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Maximum prefix length %d exceeded",sizeof(options->prefix)); 862 ierr = PetscArraycpy(options->prefix+start,prefix,n+1);CHKERRQ(ierr); 863 options->prefixstack[options->prefixind++] = start+n; 864 PetscFunctionReturn(0); 865 } 866 867 /*@C 868 PetscOptionsPrefixPop - Remove the latest options prefix, see PetscOptionsPrefixPush() for details 869 870 Logically Collective on the MPI_Comm that called PetscOptionsPrefixPush() 871 872 Input Parameters: 873 . options - options database, or NULL for the default global database 874 875 Level: advanced 876 877 .seealso: PetscOptionsPrefixPush(), PetscOptionsPush(), PetscOptionsPop(), PetscOptionsCreate(), PetscOptionsSetValue() 878 @*/ 879 PetscErrorCode PetscOptionsPrefixPop(PetscOptions options) 880 { 881 PetscInt offset; 882 883 PetscFunctionBegin; 884 options = options ? options : defaultoptions; 885 if (options->prefixind < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"More prefixes popped than pushed"); 886 options->prefixind--; 887 offset = options->prefixind ? options->prefixstack[options->prefixind-1] : 0; 888 options->prefix[offset] = 0; 889 PetscFunctionReturn(0); 890 } 891 892 /*@C 893 PetscOptionsClear - Removes all options form the database leaving it empty. 894 895 Logically Collective 896 897 Input Parameters: 898 . options - options database, use NULL for the default global database 899 900 The collectivity of this routine is complex; only the MPI processes that call this routine will 901 have the affect of these options. If some processes that create objects call this routine and others do 902 not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options 903 on different ranks. 904 905 Level: developer 906 907 .seealso: PetscOptionsInsert() 908 @*/ 909 PetscErrorCode PetscOptionsClear(PetscOptions options) 910 { 911 PetscInt i; 912 913 options = options ? options : defaultoptions; 914 if (!options) return 0; 915 916 for (i=0; i<options->N; i++) { 917 if (options->names[i]) free(options->names[i]); 918 if (options->values[i]) free(options->values[i]); 919 } 920 options->N = 0; 921 922 for (i=0; i<options->Naliases; i++) { 923 free(options->aliases1[i]); 924 free(options->aliases2[i]); 925 } 926 options->Naliases = 0; 927 928 /* destroy hash table */ 929 kh_destroy(HO,options->ht); 930 options->ht = NULL; 931 932 options->prefixind = 0; 933 options->prefix[0] = 0; 934 options->help = PETSC_FALSE; 935 return 0; 936 } 937 938 /*@C 939 PetscOptionsSetAlias - Makes a key and alias for another key 940 941 Logically Collective 942 943 Input Parameters: 944 + options - options database, or NULL for default global database 945 . newname - the alias 946 - oldname - the name that alias will refer to 947 948 Level: advanced 949 950 The collectivity of this routine is complex; only the MPI processes that call this routine will 951 have the affect of these options. If some processes that create objects call this routine and others do 952 not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options 953 on different ranks. 954 955 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(), 956 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(), 957 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 958 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 959 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 960 PetscOptionsFList(), PetscOptionsEList() 961 @*/ 962 PetscErrorCode PetscOptionsSetAlias(PetscOptions options,const char newname[],const char oldname[]) 963 { 964 PetscInt n; 965 size_t len; 966 PetscErrorCode ierr; 967 968 PetscFunctionBegin; 969 PetscValidCharPointer(newname,2); 970 PetscValidCharPointer(oldname,3); 971 options = options ? options : defaultoptions; 972 if (newname[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aliased must start with '-': Instead %s",newname); 973 if (oldname[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aliasee must start with '-': Instead %s",oldname); 974 975 n = options->Naliases; 976 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); 977 978 newname++; oldname++; 979 ierr = PetscStrlen(newname,&len);CHKERRQ(ierr); 980 options->aliases1[n] = (char*)malloc((len+1)*sizeof(char)); 981 ierr = PetscStrcpy(options->aliases1[n],newname);CHKERRQ(ierr); 982 ierr = PetscStrlen(oldname,&len);CHKERRQ(ierr); 983 options->aliases2[n] = (char*)malloc((len+1)*sizeof(char)); 984 ierr = PetscStrcpy(options->aliases2[n],oldname);CHKERRQ(ierr); 985 options->Naliases++; 986 PetscFunctionReturn(0); 987 } 988 989 /*@C 990 PetscOptionsSetValue - Sets an option name-value pair in the options 991 database, overriding whatever is already present. 992 993 Logically Collective 994 995 Input Parameters: 996 + options - options database, use NULL for the default global database 997 . name - name of option, this SHOULD have the - prepended 998 - value - the option value (not used for all options, so can be NULL) 999 1000 Level: intermediate 1001 1002 Note: 1003 This function can be called BEFORE PetscInitialize() 1004 1005 The collectivity of this routine is complex; only the MPI processes that call this routine will 1006 have the affect of these options. If some processes that create objects call this routine and others do 1007 not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options 1008 on different ranks. 1009 1010 Developers Note: Uses malloc() directly because PETSc may not be initialized yet. 1011 1012 .seealso: PetscOptionsInsert(), PetscOptionsClearValue() 1013 @*/ 1014 PetscErrorCode PetscOptionsSetValue(PetscOptions options,const char name[],const char value[]) 1015 { 1016 size_t len; 1017 int N,n,i; 1018 char **names; 1019 char fullname[MAXOPTNAME] = ""; 1020 PetscErrorCode ierr; 1021 1022 if (!options && !defaultoptions) { 1023 ierr = PetscOptionsCreateDefault();if (ierr) return ierr; 1024 } 1025 options = options ? options : defaultoptions; 1026 1027 if (name[0] != '-') return PETSC_ERR_ARG_OUTOFRANGE; 1028 1029 /* this is so that -h and -help are equivalent (p4 does not like -help)*/ 1030 if (!strcmp(name,"-h")) name = "-help"; 1031 if (!PetscOptNameCmp(name,"-help")) options->help = PETSC_TRUE; 1032 1033 name++; /* skip starting dash */ 1034 1035 if (options->prefixind > 0) { 1036 strncpy(fullname,options->prefix,sizeof(fullname)); 1037 fullname[sizeof(fullname)-1] = 0; 1038 strncat(fullname,name,sizeof(fullname)-strlen(fullname)-1); 1039 fullname[sizeof(fullname)-1] = 0; 1040 name = fullname; 1041 } 1042 1043 /* check against aliases */ 1044 N = options->Naliases; 1045 for (i=0; i<N; i++) { 1046 int result = PetscOptNameCmp(options->aliases1[i],name); 1047 if (!result) { name = options->aliases2[i]; break; } 1048 } 1049 1050 /* slow search */ 1051 N = n = options->N; 1052 names = options->names; 1053 for (i=0; i<N; i++) { 1054 int result = PetscOptNameCmp(names[i],name); 1055 if (!result) { 1056 n = i; goto setvalue; 1057 } else if (result > 0) { 1058 n = i; break; 1059 } 1060 } 1061 if (N >= MAXOPTIONS) return PETSC_ERR_MEM; 1062 /* shift remaining values up 1 */ 1063 for (i=N; i>n; i--) { 1064 options->names[i] = options->names[i-1]; 1065 options->values[i] = options->values[i-1]; 1066 options->used[i] = options->used[i-1]; 1067 } 1068 options->names[n] = NULL; 1069 options->values[n] = NULL; 1070 options->used[n] = PETSC_FALSE; 1071 options->N++; 1072 1073 /* destroy hash table */ 1074 kh_destroy(HO,options->ht); 1075 options->ht = NULL; 1076 1077 /* set new name */ 1078 len = strlen(name); 1079 options->names[n] = (char*)malloc((len+1)*sizeof(char)); 1080 if (!options->names[n]) return PETSC_ERR_MEM; 1081 strcpy(options->names[n],name); 1082 1083 setvalue: 1084 /* set new value */ 1085 if (options->values[n]) free(options->values[n]); 1086 len = value ? strlen(value) : 0; 1087 if (len) { 1088 options->values[n] = (char*)malloc((len+1)*sizeof(char)); 1089 if (!options->values[n]) return PETSC_ERR_MEM; 1090 strcpy(options->values[n],value); 1091 } else { 1092 options->values[n] = NULL; 1093 } 1094 1095 ierr = PetscOptionsMonitor(options,name,value?value:"");if (ierr) return ierr; 1096 return 0; 1097 } 1098 1099 /*@C 1100 PetscOptionsClearValue - Clears an option name-value pair in the options 1101 database, overriding whatever is already present. 1102 1103 Logically Collective 1104 1105 Input Parameter: 1106 + options - options database, use NULL for the default global database 1107 - name - name of option, this SHOULD have the - prepended 1108 1109 Level: intermediate 1110 1111 The collectivity of this routine is complex; only the MPI processes that call this routine will 1112 have the affect of these options. If some processes that create objects call this routine and others do 1113 not the code may fail in complicated ways because the same parallel solvers may incorrectly use different options 1114 on different ranks. 1115 1116 .seealso: PetscOptionsInsert() 1117 @*/ 1118 PetscErrorCode PetscOptionsClearValue(PetscOptions options,const char name[]) 1119 { 1120 int N,n,i; 1121 char **names; 1122 PetscErrorCode ierr; 1123 1124 PetscFunctionBegin; 1125 options = options ? options : defaultoptions; 1126 if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name); 1127 1128 /* this is so that -h and -help are equivalent (p4 does not like -help)*/ 1129 if (!strcmp(name,"-h")) name = "-help"; 1130 if (!PetscOptNameCmp(name,"-help")) options->help = PETSC_FALSE; 1131 1132 name++; /* skip starting dash */ 1133 1134 /* slow search */ 1135 N = n = options->N; 1136 names = options->names; 1137 for (i=0; i<N; i++) { 1138 int result = PetscOptNameCmp(names[i],name); 1139 if (!result) { 1140 n = i; break; 1141 } else if (result > 0) { 1142 n = N; break; 1143 } 1144 } 1145 if (n == N) PetscFunctionReturn(0); /* it was not present */ 1146 1147 /* remove name and value */ 1148 if (options->names[n]) free(options->names[n]); 1149 if (options->values[n]) free(options->values[n]); 1150 /* shift remaining values down 1 */ 1151 for (i=n; i<N-1; i++) { 1152 options->names[i] = options->names[i+1]; 1153 options->values[i] = options->values[i+1]; 1154 options->used[i] = options->used[i+1]; 1155 } 1156 options->N--; 1157 1158 /* destroy hash table */ 1159 kh_destroy(HO,options->ht); 1160 options->ht = NULL; 1161 1162 ierr = PetscOptionsMonitor(options,name,NULL);CHKERRQ(ierr); 1163 PetscFunctionReturn(0); 1164 } 1165 1166 /*@C 1167 PetscOptionsFindPair - Gets an option name-value pair from the options database. 1168 1169 Not Collective 1170 1171 Input Parameters: 1172 + options - options database, use NULL for the default global database 1173 . pre - the string to prepend to the name or NULL, this SHOULD NOT have the "-" prepended 1174 - name - name of option, this SHOULD have the "-" prepended 1175 1176 Output Parameters: 1177 + value - the option value (optional, not used for all options) 1178 - set - whether the option is set (optional) 1179 1180 Notes: 1181 Each process may find different values or no value depending on how options were inserted into the database 1182 1183 Level: developer 1184 1185 .seealso: PetscOptionsSetValue(), PetscOptionsClearValue() 1186 @*/ 1187 PetscErrorCode PetscOptionsFindPair(PetscOptions options,const char pre[],const char name[],const char *value[],PetscBool *set) 1188 { 1189 char buf[MAXOPTNAME]; 1190 PetscBool usehashtable = PETSC_TRUE; 1191 PetscBool matchnumbers = PETSC_TRUE; 1192 PetscErrorCode ierr; 1193 1194 PetscFunctionBegin; 1195 options = options ? options : defaultoptions; 1196 if (pre && PetscUnlikely(pre[0] == '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre); 1197 if (PetscUnlikely(name[0] != '-')) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name); 1198 1199 name++; /* skip starting dash */ 1200 1201 /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */ 1202 if (pre && pre[0]) { 1203 char *ptr = buf; 1204 if (name[0] == '-') { *ptr++ = '-'; name++; } 1205 ierr = PetscStrncpy(ptr,pre,buf+sizeof(buf)-ptr);CHKERRQ(ierr); 1206 ierr = PetscStrlcat(buf,name,sizeof(buf));CHKERRQ(ierr); 1207 name = buf; 1208 } 1209 1210 if (PetscDefined(USE_DEBUG)) { 1211 PetscBool valid; 1212 char key[MAXOPTNAME+1] = "-"; 1213 ierr = PetscStrncpy(key+1,name,sizeof(key)-1);CHKERRQ(ierr); 1214 ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr); 1215 if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name); 1216 } 1217 1218 if (!options->ht && usehashtable) { 1219 int i,ret; 1220 khiter_t it; 1221 khash_t(HO) *ht; 1222 ht = kh_init(HO); 1223 if (PetscUnlikely(!ht)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed"); 1224 ret = kh_resize(HO,ht,options->N*2); /* twice the required size to reduce risk of collisions */ 1225 if (PetscUnlikely(ret)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed"); 1226 for (i=0; i<options->N; i++) { 1227 it = kh_put(HO,ht,options->names[i],&ret); 1228 if (PetscUnlikely(ret != 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Hash table allocation failed"); 1229 kh_val(ht,it) = i; 1230 } 1231 options->ht = ht; 1232 } 1233 1234 if (usehashtable) 1235 { /* fast search */ 1236 khash_t(HO) *ht = options->ht; 1237 khiter_t it = kh_get(HO,ht,name); 1238 if (it != kh_end(ht)) { 1239 int i = kh_val(ht,it); 1240 options->used[i] = PETSC_TRUE; 1241 if (value) *value = options->values[i]; 1242 if (set) *set = PETSC_TRUE; 1243 PetscFunctionReturn(0); 1244 } 1245 } else 1246 { /* slow search */ 1247 int i, N = options->N; 1248 for (i=0; i<N; i++) { 1249 int result = PetscOptNameCmp(options->names[i],name); 1250 if (!result) { 1251 options->used[i] = PETSC_TRUE; 1252 if (value) *value = options->values[i]; 1253 if (set) *set = PETSC_TRUE; 1254 PetscFunctionReturn(0); 1255 } else if (result > 0) { 1256 break; 1257 } 1258 } 1259 } 1260 1261 /* 1262 The following block slows down all lookups in the most frequent path (most lookups are unsuccessful). 1263 Maybe this special lookup mode should be enabled on request with a push/pop API. 1264 The feature of matching _%d_ used sparingly in the codebase. 1265 */ 1266 if (matchnumbers) { 1267 int i,j,cnt = 0,locs[16],loce[16]; 1268 /* determine the location and number of all _%d_ in the key */ 1269 for (i=0; name[i]; i++) { 1270 if (name[i] == '_') { 1271 for (j=i+1; name[j]; j++) { 1272 if (name[j] >= '0' && name[j] <= '9') continue; 1273 if (name[j] == '_' && j > i+1) { /* found a number */ 1274 locs[cnt] = i+1; 1275 loce[cnt++] = j+1; 1276 } 1277 i = j-1; 1278 break; 1279 } 1280 } 1281 } 1282 for (i=0; i<cnt; i++) { 1283 PetscBool found; 1284 char opt[MAXOPTNAME+1] = "-", tmp[MAXOPTNAME]; 1285 ierr = PetscStrncpy(tmp,name,PetscMin((size_t)(locs[i]+1),sizeof(tmp)));CHKERRQ(ierr); 1286 ierr = PetscStrlcat(opt,tmp,sizeof(opt));CHKERRQ(ierr); 1287 ierr = PetscStrlcat(opt,name+loce[i],sizeof(opt));CHKERRQ(ierr); 1288 ierr = PetscOptionsFindPair(options,NULL,opt,value,&found);CHKERRQ(ierr); 1289 if (found) {if (set) *set = PETSC_TRUE; PetscFunctionReturn(0);} 1290 } 1291 } 1292 1293 if (set) *set = PETSC_FALSE; 1294 PetscFunctionReturn(0); 1295 } 1296 1297 /* Check whether any option begins with pre+name */ 1298 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions options,const char pre[], const char name[],const char *value[],PetscBool *set) 1299 { 1300 char buf[MAXOPTNAME]; 1301 int numCnt = 0, locs[16],loce[16]; 1302 PetscErrorCode ierr; 1303 1304 PetscFunctionBegin; 1305 options = options ? options : defaultoptions; 1306 if (pre && pre[0] == '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Prefix cannot begin with '-': Instead %s",pre); 1307 if (name[0] != '-') SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Name must begin with '-': Instead %s",name); 1308 1309 name++; /* skip starting dash */ 1310 1311 /* append prefix to name, if prefix="foo_" and option='--bar", prefixed option is --foo_bar */ 1312 if (pre && pre[0]) { 1313 char *ptr = buf; 1314 if (name[0] == '-') { *ptr++ = '-'; name++; } 1315 ierr = PetscStrncpy(ptr,pre,sizeof(buf)+(size_t)(ptr-buf));CHKERRQ(ierr); 1316 ierr = PetscStrlcat(buf,name,sizeof(buf));CHKERRQ(ierr); 1317 name = buf; 1318 } 1319 1320 if (PetscDefined(USE_DEBUG)) { 1321 PetscBool valid; 1322 char key[MAXOPTNAME+1] = "-"; 1323 ierr = PetscStrncpy(key+1,name,sizeof(key)-1);CHKERRQ(ierr); 1324 ierr = PetscOptionsValidKey(key,&valid);CHKERRQ(ierr); 1325 if (!valid) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid option '%s' obtained from pre='%s' and name='%s'",key,pre?pre:"",name); 1326 } 1327 1328 /* determine the location and number of all _%d_ in the key */ 1329 { 1330 int i,j; 1331 for (i=0; name[i]; i++) { 1332 if (name[i] == '_') { 1333 for (j=i+1; name[j]; j++) { 1334 if (name[j] >= '0' && name[j] <= '9') continue; 1335 if (name[j] == '_' && j > i+1) { /* found a number */ 1336 locs[numCnt] = i+1; 1337 loce[numCnt++] = j+1; 1338 } 1339 i = j-1; 1340 break; 1341 } 1342 } 1343 } 1344 } 1345 1346 { /* slow search */ 1347 int c, i; 1348 size_t len; 1349 PetscBool match; 1350 1351 for (c = -1; c < numCnt; ++c) { 1352 char opt[MAXOPTNAME+1] = "", tmp[MAXOPTNAME]; 1353 1354 if (c < 0) { 1355 ierr = PetscStrcpy(opt,name);CHKERRQ(ierr); 1356 } else { 1357 ierr = PetscStrncpy(tmp,name,PetscMin((size_t)(locs[c]+1),sizeof(tmp)));CHKERRQ(ierr); 1358 ierr = PetscStrlcat(opt,tmp,sizeof(opt));CHKERRQ(ierr); 1359 ierr = PetscStrlcat(opt,name+loce[c],sizeof(opt));CHKERRQ(ierr); 1360 } 1361 ierr = PetscStrlen(opt,&len);CHKERRQ(ierr); 1362 for (i=0; i<options->N; i++) { 1363 ierr = PetscStrncmp(options->names[i],opt,len,&match);CHKERRQ(ierr); 1364 if (match) { 1365 options->used[i] = PETSC_TRUE; 1366 if (value) *value = options->values[i]; 1367 if (set) *set = PETSC_TRUE; 1368 PetscFunctionReturn(0); 1369 } 1370 } 1371 } 1372 } 1373 1374 if (set) *set = PETSC_FALSE; 1375 PetscFunctionReturn(0); 1376 } 1377 1378 /*@C 1379 PetscOptionsReject - Generates an error if a certain option is given. 1380 1381 Not Collective 1382 1383 Input Parameters: 1384 + options - options database, use NULL for default global database 1385 . pre - the option prefix (may be NULL) 1386 . name - the option name one is seeking 1387 - mess - error message (may be NULL) 1388 1389 Level: advanced 1390 1391 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),OptionsHasName(), 1392 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1393 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1394 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1395 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1396 PetscOptionsFList(), PetscOptionsEList() 1397 @*/ 1398 PetscErrorCode PetscOptionsReject(PetscOptions options,const char pre[],const char name[],const char mess[]) 1399 { 1400 PetscErrorCode ierr; 1401 PetscBool flag = PETSC_FALSE; 1402 1403 PetscFunctionBegin; 1404 ierr = PetscOptionsHasName(options,pre,name,&flag);CHKERRQ(ierr); 1405 if (flag) { 1406 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); 1407 else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Program has disabled option: -%s%s",pre?pre:"",name+1); 1408 } 1409 PetscFunctionReturn(0); 1410 } 1411 1412 /*@C 1413 PetscOptionsHasHelp - Determines whether the "-help" option is in the database. 1414 1415 Not Collective 1416 1417 Input Parameters: 1418 . options - options database, use NULL for default global database 1419 1420 Output Parameters: 1421 . set - PETSC_TRUE if found else PETSC_FALSE. 1422 1423 Level: advanced 1424 1425 .seealso: PetscOptionsHasName() 1426 @*/ 1427 PetscErrorCode PetscOptionsHasHelp(PetscOptions options,PetscBool *set) 1428 { 1429 PetscFunctionBegin; 1430 PetscValidPointer(set,2); 1431 options = options ? options : defaultoptions; 1432 *set = options->help; 1433 PetscFunctionReturn(0); 1434 } 1435 1436 /*@C 1437 PetscOptionsHasName - Determines whether a certain option is given in the database. This returns true whether the option is a number, string or boolean, even 1438 its value is set to false. 1439 1440 Not Collective 1441 1442 Input Parameters: 1443 + options - options database, use NULL for default global database 1444 . pre - string to prepend to the name or NULL 1445 - name - the option one is seeking 1446 1447 Output Parameters: 1448 . set - PETSC_TRUE if found else PETSC_FALSE. 1449 1450 Level: beginner 1451 1452 Notes: 1453 Name cannot be simply "-h". 1454 1455 In many cases you probably want to use PetscOptionsGetBool() instead of calling this, to allowing toggling values. 1456 1457 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1458 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1459 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1460 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1461 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1462 PetscOptionsFList(), PetscOptionsEList() 1463 @*/ 1464 PetscErrorCode PetscOptionsHasName(PetscOptions options,const char pre[],const char name[],PetscBool *set) 1465 { 1466 const char *value; 1467 PetscErrorCode ierr; 1468 PetscBool flag; 1469 1470 PetscFunctionBegin; 1471 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 1472 if (set) *set = flag; 1473 PetscFunctionReturn(0); 1474 } 1475 1476 /*@C 1477 PetscOptionsGetAll - Lists all the options the program was run with in a single string. 1478 1479 Not Collective 1480 1481 Input Parameter: 1482 . options - the options database, use NULL for the default global database 1483 1484 Output Parameter: 1485 . copts - pointer where string pointer is stored 1486 1487 Notes: 1488 The array and each entry in the array should be freed with PetscFree() 1489 Each process may have different values depending on how the options were inserted into the database 1490 1491 Level: advanced 1492 1493 .seealso: PetscOptionsAllUsed(), PetscOptionsView(), PetscOptionsPush(), PetscOptionsPop() 1494 @*/ 1495 PetscErrorCode PetscOptionsGetAll(PetscOptions options,char *copts[]) 1496 { 1497 PetscErrorCode ierr; 1498 PetscInt i; 1499 size_t len = 1,lent = 0; 1500 char *coptions = NULL; 1501 1502 PetscFunctionBegin; 1503 PetscValidPointer(copts,2); 1504 options = options ? options : defaultoptions; 1505 /* count the length of the required string */ 1506 for (i=0; i<options->N; i++) { 1507 ierr = PetscStrlen(options->names[i],&lent);CHKERRQ(ierr); 1508 len += 2 + lent; 1509 if (options->values[i]) { 1510 ierr = PetscStrlen(options->values[i],&lent);CHKERRQ(ierr); 1511 len += 1 + lent; 1512 } 1513 } 1514 ierr = PetscMalloc1(len,&coptions);CHKERRQ(ierr); 1515 coptions[0] = 0; 1516 for (i=0; i<options->N; i++) { 1517 ierr = PetscStrcat(coptions,"-");CHKERRQ(ierr); 1518 ierr = PetscStrcat(coptions,options->names[i]);CHKERRQ(ierr); 1519 ierr = PetscStrcat(coptions," ");CHKERRQ(ierr); 1520 if (options->values[i]) { 1521 ierr = PetscStrcat(coptions,options->values[i]);CHKERRQ(ierr); 1522 ierr = PetscStrcat(coptions," ");CHKERRQ(ierr); 1523 } 1524 } 1525 *copts = coptions; 1526 PetscFunctionReturn(0); 1527 } 1528 1529 /*@C 1530 PetscOptionsUsed - Indicates if PETSc has used a particular option set in the database 1531 1532 Not Collective 1533 1534 Input Parameter: 1535 + options - options database, use NULL for default global database 1536 - name - string name of option 1537 1538 Output Parameter: 1539 . used - PETSC_TRUE if the option was used, otherwise false, including if option was not found in options database 1540 1541 Level: advanced 1542 1543 Notes: 1544 The value returned may be different on each process and depends on which options have been processed 1545 on the given process 1546 1547 .seealso: PetscOptionsView(), PetscOptionsLeft(), PetscOptionsAllUsed() 1548 @*/ 1549 PetscErrorCode PetscOptionsUsed(PetscOptions options,const char *name,PetscBool *used) 1550 { 1551 PetscInt i; 1552 PetscErrorCode ierr; 1553 1554 PetscFunctionBegin; 1555 PetscValidCharPointer(name,2); 1556 PetscValidPointer(used,3); 1557 options = options ? options : defaultoptions; 1558 *used = PETSC_FALSE; 1559 for (i=0; i<options->N; i++) { 1560 ierr = PetscStrcmp(options->names[i],name,used);CHKERRQ(ierr); 1561 if (*used) { 1562 *used = options->used[i]; 1563 break; 1564 } 1565 } 1566 PetscFunctionReturn(0); 1567 } 1568 1569 /*@ 1570 PetscOptionsAllUsed - Returns a count of the number of options in the 1571 database that have never been selected. 1572 1573 Not Collective 1574 1575 Input Parameter: 1576 . options - options database, use NULL for default global database 1577 1578 Output Parameter: 1579 . N - count of options not used 1580 1581 Level: advanced 1582 1583 Notes: 1584 The value returned may be different on each process and depends on which options have been processed 1585 on the given process 1586 1587 .seealso: PetscOptionsView() 1588 @*/ 1589 PetscErrorCode PetscOptionsAllUsed(PetscOptions options,PetscInt *N) 1590 { 1591 PetscInt i,n = 0; 1592 1593 PetscFunctionBegin; 1594 PetscValidIntPointer(N,2); 1595 options = options ? options : defaultoptions; 1596 for (i=0; i<options->N; i++) { 1597 if (!options->used[i]) n++; 1598 } 1599 *N = n; 1600 PetscFunctionReturn(0); 1601 } 1602 1603 /*@ 1604 PetscOptionsLeft - Prints to screen any options that were set and never used. 1605 1606 Not Collective 1607 1608 Input Parameter: 1609 . options - options database; use NULL for default global database 1610 1611 Options Database Key: 1612 . -options_left - activates PetscOptionsAllUsed() within PetscFinalize() 1613 1614 Notes: 1615 This is rarely used directly, it is called by PetscFinalize() in debug more or if -options_left 1616 is passed otherwise to help users determine possible mistakes in their usage of options. This 1617 only prints values on process zero of PETSC_COMM_WORLD. Other processes depending the objects 1618 used may have different options that are left unused. 1619 1620 Level: advanced 1621 1622 .seealso: PetscOptionsAllUsed() 1623 @*/ 1624 PetscErrorCode PetscOptionsLeft(PetscOptions options) 1625 { 1626 PetscErrorCode ierr; 1627 PetscInt i; 1628 PetscInt cnt = 0; 1629 PetscOptions toptions; 1630 1631 PetscFunctionBegin; 1632 toptions = options ? options : defaultoptions; 1633 for (i=0; i<toptions->N; i++) { 1634 if (!toptions->used[i]) { 1635 if (toptions->values[i]) { 1636 ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s value: %s\n",toptions->names[i],toptions->values[i]);CHKERRQ(ierr); 1637 } else { 1638 ierr = PetscPrintf(PETSC_COMM_WORLD,"Option left: name:-%s (no value)\n",toptions->names[i]);CHKERRQ(ierr); 1639 } 1640 } 1641 } 1642 if (!options) { 1643 toptions = defaultoptions; 1644 while (toptions->previous) { 1645 cnt++; 1646 toptions = toptions->previous; 1647 } 1648 if (cnt) { 1649 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); 1650 } 1651 } 1652 PetscFunctionReturn(0); 1653 } 1654 1655 /*@C 1656 PetscOptionsLeftGet - Returns all options that were set and never used. 1657 1658 Not Collective 1659 1660 Input Parameter: 1661 . options - options database, use NULL for default global database 1662 1663 Output Parameter: 1664 + N - count of options not used 1665 . names - names of options not used 1666 - values - values of options not used 1667 1668 Level: advanced 1669 1670 Notes: 1671 Users should call PetscOptionsLeftRestore() to free the memory allocated in this routine 1672 Notes: The value returned may be different on each process and depends on which options have been processed 1673 on the given process 1674 1675 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft() 1676 @*/ 1677 PetscErrorCode PetscOptionsLeftGet(PetscOptions options,PetscInt *N,char **names[],char **values[]) 1678 { 1679 PetscErrorCode ierr; 1680 PetscInt i,n; 1681 1682 PetscFunctionBegin; 1683 if (N) PetscValidIntPointer(N,2); 1684 if (names) PetscValidPointer(names,3); 1685 if (values) PetscValidPointer(values,4); 1686 options = options ? options : defaultoptions; 1687 1688 /* The number of unused PETSc options */ 1689 n = 0; 1690 for (i=0; i<options->N; i++) { 1691 if (!options->used[i]) n++; 1692 } 1693 if (N) { *N = n; } 1694 if (names) { ierr = PetscMalloc1(n,names);CHKERRQ(ierr); } 1695 if (values) { ierr = PetscMalloc1(n,values);CHKERRQ(ierr); } 1696 1697 n = 0; 1698 if (names || values) { 1699 for (i=0; i<options->N; i++) { 1700 if (!options->used[i]) { 1701 if (names) (*names)[n] = options->names[i]; 1702 if (values) (*values)[n] = options->values[i]; 1703 n++; 1704 } 1705 } 1706 } 1707 PetscFunctionReturn(0); 1708 } 1709 1710 /*@C 1711 PetscOptionsLeftRestore - Free memory for the unused PETSc options obtained using PetscOptionsLeftGet. 1712 1713 Not Collective 1714 1715 Input Parameter: 1716 + options - options database, use NULL for default global database 1717 . names - names of options not used 1718 - values - values of options not used 1719 1720 Level: advanced 1721 1722 .seealso: PetscOptionsAllUsed(), PetscOptionsLeft(), PetscOptionsLeftGet() 1723 @*/ 1724 PetscErrorCode PetscOptionsLeftRestore(PetscOptions options,PetscInt *N,char **names[],char **values[]) 1725 { 1726 PetscErrorCode ierr; 1727 1728 PetscFunctionBegin; 1729 if (N) PetscValidIntPointer(N,2); 1730 if (names) PetscValidPointer(names,3); 1731 if (values) PetscValidPointer(values,4); 1732 if (N) { *N = 0; } 1733 if (names) { ierr = PetscFree(*names);CHKERRQ(ierr); } 1734 if (values) { ierr = PetscFree(*values);CHKERRQ(ierr); } 1735 PetscFunctionReturn(0); 1736 } 1737 1738 /*@C 1739 PetscOptionsSetFromOptions - Sets options related to the handling of options in PETSc 1740 1741 Collective on PETSC_COMM_WORLD 1742 1743 Input Parameter: 1744 . options - options database, use NULL for default global database 1745 1746 Options Database Keys: 1747 + -options_monitor <optional filename> - prints the names and values of all runtime options as they are set. The monitor functionality is not 1748 available for options set through a file, environment variable, or on 1749 the command line. Only options set after PetscInitialize() completes will 1750 be monitored. 1751 - -options_monitor_cancel - cancel all options database monitors 1752 1753 Notes: 1754 To see all options, run your program with the -help option 1755 1756 Level: intermediate 1757 1758 @*/ 1759 PetscErrorCode PetscOptionsSetFromOptions(PetscOptions options) 1760 { 1761 PetscBool flgc = PETSC_FALSE,flgm; 1762 PetscErrorCode ierr; 1763 char monfilename[PETSC_MAX_PATH_LEN]; 1764 PetscViewer monviewer; 1765 1766 PetscFunctionBegin; 1767 /* 1768 The options argument is currently ignored since we currently maintain only a single options database 1769 1770 options = options ? options : defaultoptions; 1771 */ 1772 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for handling options","PetscOptions");CHKERRQ(ierr); 1773 ierr = PetscOptionsString("-options_monitor","Monitor options database","PetscOptionsMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flgm);CHKERRQ(ierr); 1774 ierr = PetscOptionsBool("-options_monitor_cancel","Cancel all options database monitors","PetscOptionsMonitorCancel",flgc,&flgc,NULL);CHKERRQ(ierr); 1775 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1776 if (flgm) { 1777 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,monfilename,&monviewer);CHKERRQ(ierr); 1778 ierr = PetscOptionsMonitorSet(PetscOptionsMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 1779 } 1780 if (flgc) { ierr = PetscOptionsMonitorCancel();CHKERRQ(ierr); } 1781 PetscFunctionReturn(0); 1782 } 1783 1784 /*@C 1785 PetscOptionsMonitorDefault - Print all options set value events. 1786 1787 Logically Collective on ctx 1788 1789 Input Parameters: 1790 + name - option name string 1791 . value - option value string 1792 - ctx - an ASCII viewer 1793 1794 Level: intermediate 1795 1796 Notes: 1797 The first MPI rank in the PetscViewer viewer actually prints the values, other 1798 processes may have different values set 1799 1800 .seealso: PetscOptionsMonitorSet() 1801 @*/ 1802 PetscErrorCode PetscOptionsMonitorDefault(const char name[],const char value[],void *ctx) 1803 { 1804 PetscErrorCode ierr; 1805 PetscViewer viewer = (PetscViewer)ctx; 1806 1807 PetscFunctionBegin; 1808 if (!value) { 1809 ierr = PetscViewerASCIIPrintf(viewer,"Removing option: %s\n",name,value);CHKERRQ(ierr); 1810 } else if (!value[0]) { 1811 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s (no value)\n",name);CHKERRQ(ierr); 1812 } else { 1813 ierr = PetscViewerASCIIPrintf(viewer,"Setting option: %s = %s\n",name,value);CHKERRQ(ierr); 1814 } 1815 PetscFunctionReturn(0); 1816 } 1817 1818 /*@C 1819 PetscOptionsMonitorSet - Sets an ADDITIONAL function to be called at every method that 1820 modified the PETSc options database. 1821 1822 Not Collective 1823 1824 Input Parameters: 1825 + monitor - pointer to function (if this is NULL, it turns off monitoring 1826 . mctx - [optional] context for private data for the 1827 monitor routine (use NULL if no context is desired) 1828 - monitordestroy - [optional] routine that frees monitor context 1829 (may be NULL) 1830 1831 Calling Sequence of monitor: 1832 $ monitor (const char name[], const char value[], void *mctx) 1833 1834 + name - option name string 1835 . value - option value string 1836 - mctx - optional monitoring context, as set by PetscOptionsMonitorSet() 1837 1838 Options Database Keys: 1839 + -options_monitor - sets PetscOptionsMonitorDefault() 1840 - -options_monitor_cancel - cancels all monitors that have 1841 been hardwired into a code by 1842 calls to PetscOptionsMonitorSet(), but 1843 does not cancel those set via 1844 the options database. 1845 1846 Notes: 1847 The default is to do nothing. To print the name and value of options 1848 being inserted into the database, use PetscOptionsMonitorDefault() as the monitoring routine, 1849 with a null monitoring context. 1850 1851 Several different monitoring routines may be set by calling 1852 PetscOptionsMonitorSet() multiple times; all will be called in the 1853 order in which they were set. 1854 1855 Level: beginner 1856 1857 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorCancel() 1858 @*/ 1859 PetscErrorCode PetscOptionsMonitorSet(PetscErrorCode (*monitor)(const char name[], const char value[], void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 1860 { 1861 PetscOptions options = defaultoptions; 1862 1863 PetscFunctionBegin; 1864 if (options->numbermonitors >= MAXOPTIONSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptions monitors set"); 1865 options->monitor[options->numbermonitors] = monitor; 1866 options->monitordestroy[options->numbermonitors] = monitordestroy; 1867 options->monitorcontext[options->numbermonitors++] = (void*)mctx; 1868 PetscFunctionReturn(0); 1869 } 1870 1871 /*@ 1872 PetscOptionsMonitorCancel - Clears all monitors for a PetscOptions object. 1873 1874 Not Collective 1875 1876 Options Database Key: 1877 . -options_monitor_cancel - Cancels all monitors that have 1878 been hardwired into a code by calls to PetscOptionsMonitorSet(), 1879 but does not cancel those set via the options database. 1880 1881 Level: intermediate 1882 1883 .seealso: PetscOptionsMonitorDefault(), PetscOptionsMonitorSet() 1884 @*/ 1885 PetscErrorCode PetscOptionsMonitorCancel(void) 1886 { 1887 PetscErrorCode ierr; 1888 PetscInt i; 1889 PetscOptions options = defaultoptions; 1890 1891 PetscFunctionBegin; 1892 for (i=0; i<options->numbermonitors; i++) { 1893 if (options->monitordestroy[i]) { 1894 ierr = (*options->monitordestroy[i])(&options->monitorcontext[i]);CHKERRQ(ierr); 1895 } 1896 } 1897 options->numbermonitors = 0; 1898 PetscFunctionReturn(0); 1899 } 1900 1901 /* 1902 PetscOptionsStringToBool - Converts string to PetscBool , handles cases like "yes", "no", "true", "false", "0", "1", "off", "on". 1903 */ 1904 PetscErrorCode PetscOptionsStringToBool(const char value[],PetscBool *a) 1905 { 1906 PetscBool istrue,isfalse; 1907 size_t len; 1908 PetscErrorCode ierr; 1909 1910 PetscFunctionBegin; 1911 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1912 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Character string of length zero has no logical value"); 1913 ierr = PetscStrcasecmp(value,"TRUE",&istrue);CHKERRQ(ierr); 1914 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1915 ierr = PetscStrcasecmp(value,"YES",&istrue);CHKERRQ(ierr); 1916 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1917 ierr = PetscStrcasecmp(value,"1",&istrue);CHKERRQ(ierr); 1918 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1919 ierr = PetscStrcasecmp(value,"on",&istrue);CHKERRQ(ierr); 1920 if (istrue) {*a = PETSC_TRUE; PetscFunctionReturn(0);} 1921 ierr = PetscStrcasecmp(value,"FALSE",&isfalse);CHKERRQ(ierr); 1922 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1923 ierr = PetscStrcasecmp(value,"NO",&isfalse);CHKERRQ(ierr); 1924 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1925 ierr = PetscStrcasecmp(value,"0",&isfalse);CHKERRQ(ierr); 1926 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1927 ierr = PetscStrcasecmp(value,"off",&isfalse);CHKERRQ(ierr); 1928 if (isfalse) {*a = PETSC_FALSE; PetscFunctionReturn(0);} 1929 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown logical value: %s",value); 1930 } 1931 1932 /* 1933 PetscOptionsStringToInt - Converts a string to an integer value. Handles special cases such as "default" and "decide" 1934 */ 1935 PetscErrorCode PetscOptionsStringToInt(const char name[],PetscInt *a) 1936 { 1937 PetscErrorCode ierr; 1938 size_t len; 1939 PetscBool decide,tdefault,mouse; 1940 1941 PetscFunctionBegin; 1942 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 1943 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 1944 1945 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&tdefault);CHKERRQ(ierr); 1946 if (!tdefault) { 1947 ierr = PetscStrcasecmp(name,"DEFAULT",&tdefault);CHKERRQ(ierr); 1948 } 1949 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&decide);CHKERRQ(ierr); 1950 if (!decide) { 1951 ierr = PetscStrcasecmp(name,"DECIDE",&decide);CHKERRQ(ierr); 1952 } 1953 ierr = PetscStrcasecmp(name,"mouse",&mouse);CHKERRQ(ierr); 1954 1955 if (tdefault) *a = PETSC_DEFAULT; 1956 else if (decide) *a = PETSC_DECIDE; 1957 else if (mouse) *a = -1; 1958 else { 1959 char *endptr; 1960 long strtolval; 1961 1962 strtolval = strtol(name,&endptr,10); 1963 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); 1964 1965 #if defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE_ATOLL) 1966 (void) strtolval; 1967 *a = atoll(name); 1968 #elif defined(PETSC_USE_64BIT_INDICES) && defined(PETSC_HAVE___INT64) 1969 (void) strtolval; 1970 *a = _atoi64(name); 1971 #else 1972 *a = (PetscInt)strtolval; 1973 #endif 1974 } 1975 PetscFunctionReturn(0); 1976 } 1977 1978 #if defined(PETSC_USE_REAL___FLOAT128) 1979 #include <quadmath.h> 1980 #endif 1981 1982 static PetscErrorCode PetscStrtod(const char name[],PetscReal *a,char **endptr) 1983 { 1984 PetscFunctionBegin; 1985 #if defined(PETSC_USE_REAL___FLOAT128) 1986 *a = strtoflt128(name,endptr); 1987 #else 1988 *a = (PetscReal)strtod(name,endptr); 1989 #endif 1990 PetscFunctionReturn(0); 1991 } 1992 1993 static PetscErrorCode PetscStrtoz(const char name[],PetscScalar *a,char **endptr,PetscBool *isImaginary) 1994 { 1995 PetscBool hasi = PETSC_FALSE; 1996 char *ptr; 1997 PetscReal strtoval; 1998 PetscErrorCode ierr; 1999 2000 PetscFunctionBegin; 2001 ierr = PetscStrtod(name,&strtoval,&ptr);CHKERRQ(ierr); 2002 if (ptr == name) { 2003 strtoval = 1.; 2004 hasi = PETSC_TRUE; 2005 if (name[0] == 'i') { 2006 ptr++; 2007 } else if (name[0] == '+' && name[1] == 'i') { 2008 ptr += 2; 2009 } else if (name[0] == '-' && name[1] == 'i') { 2010 strtoval = -1.; 2011 ptr += 2; 2012 } 2013 } else if (*ptr == 'i') { 2014 hasi = PETSC_TRUE; 2015 ptr++; 2016 } 2017 *endptr = ptr; 2018 *isImaginary = hasi; 2019 if (hasi) { 2020 #if !defined(PETSC_USE_COMPLEX) 2021 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s contains imaginary but complex not supported ",name); 2022 #else 2023 *a = PetscCMPLX(0.,strtoval); 2024 #endif 2025 } else { 2026 *a = strtoval; 2027 } 2028 PetscFunctionReturn(0); 2029 } 2030 2031 /* 2032 Converts a string to PetscReal value. Handles special cases like "default" and "decide" 2033 */ 2034 PetscErrorCode PetscOptionsStringToReal(const char name[],PetscReal *a) 2035 { 2036 size_t len; 2037 PetscBool match; 2038 char *endptr; 2039 PetscErrorCode ierr; 2040 2041 PetscFunctionBegin; 2042 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 2043 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"String of length zero has no numerical value"); 2044 2045 ierr = PetscStrcasecmp(name,"PETSC_DEFAULT",&match);CHKERRQ(ierr); 2046 if (!match) { 2047 ierr = PetscStrcasecmp(name,"DEFAULT",&match);CHKERRQ(ierr); 2048 } 2049 if (match) {*a = PETSC_DEFAULT; PetscFunctionReturn(0);} 2050 2051 ierr = PetscStrcasecmp(name,"PETSC_DECIDE",&match);CHKERRQ(ierr); 2052 if (!match) { 2053 ierr = PetscStrcasecmp(name,"DECIDE",&match);CHKERRQ(ierr); 2054 } 2055 if (match) {*a = PETSC_DECIDE; PetscFunctionReturn(0);} 2056 2057 ierr = PetscStrtod(name,a,&endptr);CHKERRQ(ierr); 2058 if ((size_t) (endptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value",name); 2059 PetscFunctionReturn(0); 2060 } 2061 2062 PetscErrorCode PetscOptionsStringToScalar(const char name[],PetscScalar *a) 2063 { 2064 PetscBool imag1; 2065 size_t len; 2066 PetscScalar val = 0.; 2067 char *ptr = NULL; 2068 PetscErrorCode ierr; 2069 2070 PetscFunctionBegin; 2071 ierr = PetscStrlen(name,&len);CHKERRQ(ierr); 2072 if (!len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"character string of length zero has no numerical value"); 2073 ierr = PetscStrtoz(name,&val,&ptr,&imag1);CHKERRQ(ierr); 2074 #if defined(PETSC_USE_COMPLEX) 2075 if ((size_t) (ptr - name) < len) { 2076 PetscBool imag2; 2077 PetscScalar val2; 2078 2079 ierr = PetscStrtoz(ptr,&val2,&ptr,&imag2);CHKERRQ(ierr); 2080 if (imag1 || !imag2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s: must specify imaginary component second",name); 2081 val = PetscCMPLX(PetscRealPart(val),PetscImaginaryPart(val2)); 2082 } 2083 #endif 2084 if ((size_t) (ptr - name) != len) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Input string %s has no numeric value ",name); 2085 *a = val; 2086 PetscFunctionReturn(0); 2087 } 2088 2089 /*@C 2090 PetscOptionsGetBool - Gets the Logical (true or false) value for a particular 2091 option in the database. 2092 2093 Not Collective 2094 2095 Input Parameters: 2096 + options - options database, use NULL for default global database 2097 . pre - the string to prepend to the name or NULL 2098 - name - the option one is seeking 2099 2100 Output Parameter: 2101 + ivalue - the logical value to return 2102 - set - PETSC_TRUE if found, else PETSC_FALSE 2103 2104 Level: beginner 2105 2106 Notes: 2107 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 2108 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 2109 2110 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 2111 is equivalent to -requested_bool true 2112 2113 If the user does not supply the option at all ivalue is NOT changed. Thus 2114 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2115 2116 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 2117 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsGetInt(), PetscOptionsBool(), 2118 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2119 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2120 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2121 PetscOptionsFList(), PetscOptionsEList() 2122 @*/ 2123 PetscErrorCode PetscOptionsGetBool(PetscOptions options,const char pre[],const char name[],PetscBool *ivalue,PetscBool *set) 2124 { 2125 const char *value; 2126 PetscBool flag; 2127 PetscErrorCode ierr; 2128 2129 PetscFunctionBegin; 2130 PetscValidCharPointer(name,3); 2131 if (ivalue) PetscValidIntPointer(ivalue,4); 2132 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2133 if (flag) { 2134 if (set) *set = PETSC_TRUE; 2135 if (!value) { 2136 if (ivalue) *ivalue = PETSC_TRUE; 2137 } else { 2138 ierr = PetscOptionsStringToBool(value, &flag);CHKERRQ(ierr); 2139 if (ivalue) *ivalue = flag; 2140 } 2141 } else { 2142 if (set) *set = PETSC_FALSE; 2143 } 2144 PetscFunctionReturn(0); 2145 } 2146 2147 /*@C 2148 PetscOptionsGetEList - Puts a list of option values that a single one may be selected from 2149 2150 Not Collective 2151 2152 Input Parameters: 2153 + options - options database, use NULL for default global database 2154 . pre - the string to prepend to the name or NULL 2155 . opt - option name 2156 . list - the possible choices (one of these must be selected, anything else is invalid) 2157 - ntext - number of choices 2158 2159 Output Parameter: 2160 + value - the index of the value to return (defaults to zero if the option name is given but no choice is listed) 2161 - set - PETSC_TRUE if found, else PETSC_FALSE 2162 2163 Level: intermediate 2164 2165 Notes: 2166 If the user does not supply the option value is NOT changed. Thus 2167 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2168 2169 See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 2170 2171 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2172 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2173 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2174 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2175 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2176 PetscOptionsFList(), PetscOptionsEList() 2177 @*/ 2178 PetscErrorCode PetscOptionsGetEList(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscInt ntext,PetscInt *value,PetscBool *set) 2179 { 2180 PetscErrorCode ierr; 2181 size_t alen,len = 0, tlen = 0; 2182 char *svalue; 2183 PetscBool aset,flg = PETSC_FALSE; 2184 PetscInt i; 2185 2186 PetscFunctionBegin; 2187 PetscValidCharPointer(opt,3); 2188 for (i=0; i<ntext; i++) { 2189 ierr = PetscStrlen(list[i],&alen);CHKERRQ(ierr); 2190 if (alen > len) len = alen; 2191 tlen += len + 1; 2192 } 2193 len += 5; /* a little extra space for user mistypes */ 2194 ierr = PetscMalloc1(len,&svalue);CHKERRQ(ierr); 2195 ierr = PetscOptionsGetString(options,pre,opt,svalue,len,&aset);CHKERRQ(ierr); 2196 if (aset) { 2197 ierr = PetscEListFind(ntext,list,svalue,value,&flg);CHKERRQ(ierr); 2198 if (!flg) { 2199 char *avail,*pavl; 2200 2201 ierr = PetscMalloc1(tlen,&avail);CHKERRQ(ierr); 2202 pavl = avail; 2203 for (i=0; i<ntext; i++) { 2204 ierr = PetscStrlen(list[i],&alen);CHKERRQ(ierr); 2205 ierr = PetscStrcpy(pavl,list[i]);CHKERRQ(ierr); 2206 pavl += alen; 2207 ierr = PetscStrcpy(pavl," ");CHKERRQ(ierr); 2208 pavl += 1; 2209 } 2210 ierr = PetscStrtolower(avail);CHKERRQ(ierr); 2211 SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown option %s for -%s%s. Available options: %s",svalue,pre ? pre : "",opt+1,avail); 2212 } 2213 if (set) *set = PETSC_TRUE; 2214 } else if (set) *set = PETSC_FALSE; 2215 ierr = PetscFree(svalue);CHKERRQ(ierr); 2216 PetscFunctionReturn(0); 2217 } 2218 2219 /*@C 2220 PetscOptionsGetEnum - Gets the enum value for a particular option in the database. 2221 2222 Not Collective 2223 2224 Input Parameters: 2225 + options - options database, use NULL for default global database 2226 . pre - option prefix or NULL 2227 . opt - option name 2228 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2229 - defaultv - the default (current) value 2230 2231 Output Parameter: 2232 + value - the value to return 2233 - set - PETSC_TRUE if found, else PETSC_FALSE 2234 2235 Level: beginner 2236 2237 Notes: 2238 If the user does not supply the option value is NOT changed. Thus 2239 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2240 2241 List is usually something like PCASMTypes or some other predefined list of enum names 2242 2243 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2244 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2245 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2246 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2247 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2248 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2249 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2250 @*/ 2251 PetscErrorCode PetscOptionsGetEnum(PetscOptions options,const char pre[],const char opt[],const char * const *list,PetscEnum *value,PetscBool *set) 2252 { 2253 PetscErrorCode ierr; 2254 PetscInt ntext = 0,tval; 2255 PetscBool fset; 2256 2257 PetscFunctionBegin; 2258 PetscValidCharPointer(opt,3); 2259 while (list[ntext++]) { 2260 if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 2261 } 2262 if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 2263 ntext -= 3; 2264 ierr = PetscOptionsGetEList(options,pre,opt,list,ntext,&tval,&fset);CHKERRQ(ierr); 2265 /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 2266 if (fset) *value = (PetscEnum)tval; 2267 if (set) *set = fset; 2268 PetscFunctionReturn(0); 2269 } 2270 2271 /*@C 2272 PetscOptionsGetInt - Gets the integer value for a particular option in the database. 2273 2274 Not Collective 2275 2276 Input Parameters: 2277 + options - options database, use NULL for default global database 2278 . pre - the string to prepend to the name or NULL 2279 - name - the option one is seeking 2280 2281 Output Parameter: 2282 + ivalue - the integer value to return 2283 - set - PETSC_TRUE if found, else PETSC_FALSE 2284 2285 Level: beginner 2286 2287 Notes: 2288 If the user does not supply the option ivalue is NOT changed. Thus 2289 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2290 2291 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 2292 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2293 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 2294 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2295 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2296 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2297 PetscOptionsFList(), PetscOptionsEList() 2298 @*/ 2299 PetscErrorCode PetscOptionsGetInt(PetscOptions options,const char pre[],const char name[],PetscInt *ivalue,PetscBool *set) 2300 { 2301 const char *value; 2302 PetscErrorCode ierr; 2303 PetscBool flag; 2304 2305 PetscFunctionBegin; 2306 PetscValidCharPointer(name,3); 2307 PetscValidIntPointer(ivalue,4); 2308 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2309 if (flag) { 2310 if (!value) { 2311 if (set) *set = PETSC_FALSE; 2312 } else { 2313 if (set) *set = PETSC_TRUE; 2314 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2315 } 2316 } else { 2317 if (set) *set = PETSC_FALSE; 2318 } 2319 PetscFunctionReturn(0); 2320 } 2321 2322 /*@C 2323 PetscOptionsGetReal - Gets the double precision value for a particular 2324 option in the database. 2325 2326 Not Collective 2327 2328 Input Parameters: 2329 + options - options database, use NULL for default global database 2330 . pre - string to prepend to each name or NULL 2331 - name - the option one is seeking 2332 2333 Output Parameter: 2334 + dvalue - the double value to return 2335 - set - PETSC_TRUE if found, PETSC_FALSE if not found 2336 2337 Notes: 2338 If the user does not supply the option dvalue is NOT changed. Thus 2339 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2340 2341 Level: beginner 2342 2343 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2344 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(),PetscOptionsBool(), 2345 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2346 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2347 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2348 PetscOptionsFList(), PetscOptionsEList() 2349 @*/ 2350 PetscErrorCode PetscOptionsGetReal(PetscOptions options,const char pre[],const char name[],PetscReal *dvalue,PetscBool *set) 2351 { 2352 const char *value; 2353 PetscBool flag; 2354 PetscErrorCode ierr; 2355 2356 PetscFunctionBegin; 2357 PetscValidCharPointer(name,3); 2358 PetscValidRealPointer(dvalue,4); 2359 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2360 if (flag) { 2361 if (!value) { 2362 if (set) *set = PETSC_FALSE; 2363 } else { 2364 if (set) *set = PETSC_TRUE; 2365 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2366 } 2367 } else { 2368 if (set) *set = PETSC_FALSE; 2369 } 2370 PetscFunctionReturn(0); 2371 } 2372 2373 /*@C 2374 PetscOptionsGetScalar - Gets the scalar value for a particular 2375 option in the database. 2376 2377 Not Collective 2378 2379 Input Parameters: 2380 + options - options database, use NULL for default global database 2381 . pre - string to prepend to each name or NULL 2382 - name - the option one is seeking 2383 2384 Output Parameter: 2385 + dvalue - the double value to return 2386 - set - PETSC_TRUE if found, else PETSC_FALSE 2387 2388 Level: beginner 2389 2390 Usage: 2391 A complex number 2+3i must be specified with NO spaces 2392 2393 Notes: 2394 If the user does not supply the option dvalue is NOT changed. Thus 2395 you should ALWAYS initialize the ivalue if you access it without first checking if the set flag is true. 2396 2397 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2398 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2399 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2400 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2401 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2402 PetscOptionsFList(), PetscOptionsEList() 2403 @*/ 2404 PetscErrorCode PetscOptionsGetScalar(PetscOptions options,const char pre[],const char name[],PetscScalar *dvalue,PetscBool *set) 2405 { 2406 const char *value; 2407 PetscBool flag; 2408 PetscErrorCode ierr; 2409 2410 PetscFunctionBegin; 2411 PetscValidCharPointer(name,3); 2412 PetscValidScalarPointer(dvalue,4); 2413 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2414 if (flag) { 2415 if (!value) { 2416 if (set) *set = PETSC_FALSE; 2417 } else { 2418 #if !defined(PETSC_USE_COMPLEX) 2419 ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 2420 #else 2421 ierr = PetscOptionsStringToScalar(value,dvalue);CHKERRQ(ierr); 2422 #endif 2423 if (set) *set = PETSC_TRUE; 2424 } 2425 } else { /* flag */ 2426 if (set) *set = PETSC_FALSE; 2427 } 2428 PetscFunctionReturn(0); 2429 } 2430 2431 /*@C 2432 PetscOptionsGetString - Gets the string value for a particular option in 2433 the database. 2434 2435 Not Collective 2436 2437 Input Parameters: 2438 + options - options database, use NULL for default global database 2439 . pre - string to prepend to name or NULL 2440 . name - the option one is seeking 2441 - len - maximum length of the string including null termination 2442 2443 Output Parameters: 2444 + string - location to copy string 2445 - set - PETSC_TRUE if found, else PETSC_FALSE 2446 2447 Level: beginner 2448 2449 Fortran Note: 2450 The Fortran interface is slightly different from the C/C++ 2451 interface (len is not used). Sample usage in Fortran follows 2452 .vb 2453 character *20 string 2454 PetscErrorCode ierr 2455 PetscBool set 2456 call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s',string,set,ierr) 2457 .ve 2458 2459 Notes: 2460 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 2461 2462 If the user does not use the option then the string is not changed. Thus 2463 you should ALWAYS initialize the string if you access it without first checking if the set flag is true. 2464 2465 Note: 2466 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). 2467 2468 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2469 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2470 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2471 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2472 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2473 PetscOptionsFList(), PetscOptionsEList() 2474 @*/ 2475 PetscErrorCode PetscOptionsGetString(PetscOptions options,const char pre[],const char name[],char string[],size_t len,PetscBool *set) 2476 { 2477 const char *value; 2478 PetscBool flag; 2479 PetscErrorCode ierr; 2480 2481 PetscFunctionBegin; 2482 PetscValidCharPointer(name,3); 2483 PetscValidCharPointer(string,4); 2484 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);CHKERRQ(ierr); 2485 if (!flag) { 2486 if (set) *set = PETSC_FALSE; 2487 } else { 2488 if (set) *set = PETSC_TRUE; 2489 if (value) { 2490 ierr = PetscStrncpy(string,value,len);CHKERRQ(ierr); 2491 } else { 2492 ierr = PetscArrayzero(string,len);CHKERRQ(ierr); 2493 } 2494 } 2495 PetscFunctionReturn(0); 2496 } 2497 2498 char *PetscOptionsGetStringMatlab(PetscOptions options,const char pre[],const char name[]) 2499 { 2500 const char *value; 2501 PetscBool flag; 2502 PetscErrorCode ierr; 2503 2504 PetscFunctionBegin; 2505 ierr = PetscOptionsFindPair(options,pre,name,&value,&flag);if (ierr) PetscFunctionReturn(NULL); 2506 if (flag) PetscFunctionReturn((char*)value); 2507 else PetscFunctionReturn(NULL); 2508 } 2509 2510 /*@C 2511 PetscOptionsGetBoolArray - Gets an array of Logical (true or false) values for a particular 2512 option in the database. The values must be separated with commas with 2513 no intervening spaces. 2514 2515 Not Collective 2516 2517 Input Parameters: 2518 + options - options database, use NULL for default global database 2519 . pre - string to prepend to each name or NULL 2520 . name - the option one is seeking 2521 - nmax - maximum number of values to retrieve 2522 2523 Output Parameter: 2524 + dvalue - the integer values to return 2525 . nmax - actual number of values retreived 2526 - set - PETSC_TRUE if found, else PETSC_FALSE 2527 2528 Level: beginner 2529 2530 Notes: 2531 TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 2532 FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 2533 2534 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2535 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2536 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2537 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2538 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2539 PetscOptionsFList(), PetscOptionsEList() 2540 @*/ 2541 PetscErrorCode PetscOptionsGetBoolArray(PetscOptions options,const char pre[],const char name[],PetscBool dvalue[],PetscInt *nmax,PetscBool *set) 2542 { 2543 const char *svalue; 2544 char *value; 2545 PetscErrorCode ierr; 2546 PetscInt n = 0; 2547 PetscBool flag; 2548 PetscToken token; 2549 2550 PetscFunctionBegin; 2551 PetscValidCharPointer(name,3); 2552 PetscValidIntPointer(dvalue,4); 2553 PetscValidIntPointer(nmax,5); 2554 2555 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2556 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2557 if (set) *set = PETSC_TRUE; 2558 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2559 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2560 while (value && n < *nmax) { 2561 ierr = PetscOptionsStringToBool(value,dvalue);CHKERRQ(ierr); 2562 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2563 dvalue++; 2564 n++; 2565 } 2566 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2567 *nmax = n; 2568 PetscFunctionReturn(0); 2569 } 2570 2571 /*@C 2572 PetscOptionsGetEnumArray - Gets an array of enum values for a particular option in the database. 2573 2574 Not Collective 2575 2576 Input Parameters: 2577 + options - options database, use NULL for default global database 2578 . pre - option prefix or NULL 2579 . name - option name 2580 . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 2581 - nmax - maximum number of values to retrieve 2582 2583 Output Parameters: 2584 + ivalue - the enum values to return 2585 . nmax - actual number of values retreived 2586 - set - PETSC_TRUE if found, else PETSC_FALSE 2587 2588 Level: beginner 2589 2590 Notes: 2591 The array must be passed as a comma separated list. 2592 2593 There must be no intervening spaces between the values. 2594 2595 list is usually something like PCASMTypes or some other predefined list of enum names. 2596 2597 .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 2598 PetscOptionsGetEnum(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 2599 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), PetscOptionsName(), 2600 PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), PetscOptionsStringArray(),PetscOptionsRealArray(), 2601 PetscOptionsScalar(), PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2602 PetscOptionsFList(), PetscOptionsEList(), PetscOptionsGetEList(), PetscOptionsEnum() 2603 @*/ 2604 PetscErrorCode PetscOptionsGetEnumArray(PetscOptions options,const char pre[],const char name[],const char *const *list,PetscEnum ivalue[],PetscInt *nmax,PetscBool *set) 2605 { 2606 const char *svalue; 2607 char *value; 2608 PetscInt n = 0; 2609 PetscEnum evalue; 2610 PetscBool flag; 2611 PetscToken token; 2612 PetscErrorCode ierr; 2613 2614 PetscFunctionBegin; 2615 PetscValidCharPointer(name,3); 2616 PetscValidPointer(list,4); 2617 PetscValidPointer(ivalue,5); 2618 PetscValidIntPointer(nmax,6); 2619 2620 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2621 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2622 if (set) *set = PETSC_TRUE; 2623 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2624 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2625 while (value && n < *nmax) { 2626 ierr = PetscEnumFind(list,value,&evalue,&flag);CHKERRQ(ierr); 2627 if (!flag) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown enum value '%s' for -%s%s",svalue,pre ? pre : "",name+1); 2628 ivalue[n++] = evalue; 2629 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2630 } 2631 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2632 *nmax = n; 2633 PetscFunctionReturn(0); 2634 } 2635 2636 /*@C 2637 PetscOptionsGetIntArray - Gets an array of integer values for a particular 2638 option in the database. 2639 2640 Not Collective 2641 2642 Input Parameters: 2643 + options - options database, use NULL for default global database 2644 . pre - string to prepend to each name or NULL 2645 . name - the option one is seeking 2646 - nmax - maximum number of values to retrieve 2647 2648 Output Parameter: 2649 + ivalue - the integer values to return 2650 . nmax - actual number of values retreived 2651 - set - PETSC_TRUE if found, else PETSC_FALSE 2652 2653 Level: beginner 2654 2655 Notes: 2656 The array can be passed as 2657 a comma separated list: 0,1,2,3,4,5,6,7 2658 a range (start-end+1): 0-8 2659 a range with given increment (start-end+1:inc): 0-7:2 2660 a combination of values and ranges separated by commas: 0,1-8,8-15:2 2661 2662 There must be no intervening spaces between the values. 2663 2664 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2665 PetscOptionsGetString(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2666 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2667 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2668 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2669 PetscOptionsFList(), PetscOptionsEList() 2670 @*/ 2671 PetscErrorCode PetscOptionsGetIntArray(PetscOptions options,const char pre[],const char name[],PetscInt ivalue[],PetscInt *nmax,PetscBool *set) 2672 { 2673 const char *svalue; 2674 char *value; 2675 PetscErrorCode ierr; 2676 PetscInt n = 0,i,j,start,end,inc,nvalues; 2677 size_t len; 2678 PetscBool flag,foundrange; 2679 PetscToken token; 2680 2681 PetscFunctionBegin; 2682 PetscValidCharPointer(name,3); 2683 PetscValidIntPointer(ivalue,4); 2684 PetscValidIntPointer(nmax,5); 2685 2686 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2687 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2688 if (set) *set = PETSC_TRUE; 2689 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2690 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2691 while (value && n < *nmax) { 2692 /* look for form d-D where d and D are integers */ 2693 foundrange = PETSC_FALSE; 2694 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 2695 if (value[0] == '-') i=2; 2696 else i=1; 2697 for (;i<(int)len; i++) { 2698 if (value[i] == '-') { 2699 if (i == (int)len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 2700 value[i] = 0; 2701 2702 ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 2703 inc = 1; 2704 j = i+1; 2705 for (;j<(int)len; j++) { 2706 if (value[j] == ':') { 2707 value[j] = 0; 2708 2709 ierr = PetscOptionsStringToInt(value+j+1,&inc);CHKERRQ(ierr); 2710 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); 2711 break; 2712 } 2713 } 2714 ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 2715 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); 2716 nvalues = (end-start)/inc + (end-start)%inc; 2717 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); 2718 for (;start<end; start+=inc) { 2719 *ivalue = start; ivalue++;n++; 2720 } 2721 foundrange = PETSC_TRUE; 2722 break; 2723 } 2724 } 2725 if (!foundrange) { 2726 ierr = PetscOptionsStringToInt(value,ivalue);CHKERRQ(ierr); 2727 ivalue++; 2728 n++; 2729 } 2730 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2731 } 2732 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2733 *nmax = n; 2734 PetscFunctionReturn(0); 2735 } 2736 2737 /*@C 2738 PetscOptionsGetRealArray - Gets an array of double precision values for a 2739 particular option in the database. The values must be separated with 2740 commas with no intervening spaces. 2741 2742 Not Collective 2743 2744 Input Parameters: 2745 + options - options database, use NULL for default global database 2746 . pre - string to prepend to each name or NULL 2747 . name - the option one is seeking 2748 - nmax - maximum number of values to retrieve 2749 2750 Output Parameters: 2751 + dvalue - the double values to return 2752 . nmax - actual number of values retreived 2753 - set - PETSC_TRUE if found, else PETSC_FALSE 2754 2755 Level: beginner 2756 2757 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2758 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2759 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2760 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2761 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2762 PetscOptionsFList(), PetscOptionsEList() 2763 @*/ 2764 PetscErrorCode PetscOptionsGetRealArray(PetscOptions options,const char pre[],const char name[],PetscReal dvalue[],PetscInt *nmax,PetscBool *set) 2765 { 2766 const char *svalue; 2767 char *value; 2768 PetscErrorCode ierr; 2769 PetscInt n = 0; 2770 PetscBool flag; 2771 PetscToken token; 2772 2773 PetscFunctionBegin; 2774 PetscValidCharPointer(name,3); 2775 PetscValidRealPointer(dvalue,4); 2776 PetscValidIntPointer(nmax,5); 2777 2778 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2779 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2780 if (set) *set = PETSC_TRUE; 2781 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2782 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2783 while (value && n < *nmax) { 2784 ierr = PetscOptionsStringToReal(value,dvalue++);CHKERRQ(ierr); 2785 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2786 n++; 2787 } 2788 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2789 *nmax = n; 2790 PetscFunctionReturn(0); 2791 } 2792 2793 /*@C 2794 PetscOptionsGetScalarArray - Gets an array of scalars for a 2795 particular option in the database. The values must be separated with 2796 commas with no intervening spaces. 2797 2798 Not Collective 2799 2800 Input Parameters: 2801 + options - options database, use NULL for default global database 2802 . pre - string to prepend to each name or NULL 2803 . name - the option one is seeking 2804 - nmax - maximum number of values to retrieve 2805 2806 Output Parameters: 2807 + dvalue - the scalar values to return 2808 . nmax - actual number of values retreived 2809 - set - PETSC_TRUE if found, else PETSC_FALSE 2810 2811 Level: beginner 2812 2813 .seealso: PetscOptionsGetInt(), PetscOptionsHasName(), 2814 PetscOptionsGetString(), PetscOptionsGetIntArray(), PetscOptionsBool(), 2815 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2816 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2817 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2818 PetscOptionsFList(), PetscOptionsEList() 2819 @*/ 2820 PetscErrorCode PetscOptionsGetScalarArray(PetscOptions options,const char pre[],const char name[],PetscScalar dvalue[],PetscInt *nmax,PetscBool *set) 2821 { 2822 const char *svalue; 2823 char *value; 2824 PetscErrorCode ierr; 2825 PetscInt n = 0; 2826 PetscBool flag; 2827 PetscToken token; 2828 2829 PetscFunctionBegin; 2830 PetscValidCharPointer(name,3); 2831 PetscValidRealPointer(dvalue,4); 2832 PetscValidIntPointer(nmax,5); 2833 2834 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2835 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2836 if (set) *set = PETSC_TRUE; 2837 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2838 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2839 while (value && n < *nmax) { 2840 ierr = PetscOptionsStringToScalar(value,dvalue++);CHKERRQ(ierr); 2841 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2842 n++; 2843 } 2844 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2845 *nmax = n; 2846 PetscFunctionReturn(0); 2847 } 2848 2849 /*@C 2850 PetscOptionsGetStringArray - Gets an array of string values for a particular 2851 option in the database. The values must be separated with commas with 2852 no intervening spaces. 2853 2854 Not Collective 2855 2856 Input Parameters: 2857 + options - options database, use NULL for default global database 2858 . pre - string to prepend to name or NULL 2859 . name - the option one is seeking 2860 - nmax - maximum number of strings 2861 2862 Output Parameter: 2863 + strings - location to copy strings 2864 - set - PETSC_TRUE if found, else PETSC_FALSE 2865 2866 Level: beginner 2867 2868 Notes: 2869 The user should pass in an array of pointers to char, to hold all the 2870 strings returned by this function. 2871 2872 The user is responsible for deallocating the strings that are 2873 returned. The Fortran interface for this routine is not supported. 2874 2875 .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 2876 PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 2877 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 2878 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 2879 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 2880 PetscOptionsFList(), PetscOptionsEList() 2881 @*/ 2882 PetscErrorCode PetscOptionsGetStringArray(PetscOptions options,const char pre[],const char name[],char *strings[],PetscInt *nmax,PetscBool *set) 2883 { 2884 const char *svalue; 2885 char *value; 2886 PetscErrorCode ierr; 2887 PetscInt n = 0; 2888 PetscBool flag; 2889 PetscToken token; 2890 2891 PetscFunctionBegin; 2892 PetscValidCharPointer(name,3); 2893 PetscValidPointer(strings,4); 2894 PetscValidIntPointer(nmax,5); 2895 2896 ierr = PetscOptionsFindPair(options,pre,name,&svalue,&flag);CHKERRQ(ierr); 2897 if (!flag || !svalue) { if (set) *set = PETSC_FALSE; *nmax = 0; PetscFunctionReturn(0);} 2898 if (set) *set = PETSC_TRUE; 2899 ierr = PetscTokenCreate(svalue,',',&token);CHKERRQ(ierr); 2900 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2901 while (value && n < *nmax) { 2902 ierr = PetscStrallocpy(value,&strings[n]);CHKERRQ(ierr); 2903 ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 2904 n++; 2905 } 2906 ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 2907 *nmax = n; 2908 PetscFunctionReturn(0); 2909 } 2910 2911 /*@C 2912 PetscOptionsDeprecated - mark an option as deprecated, optionally replacing it with a new one 2913 2914 Prints a deprecation warning, unless an option is supplied to suppress. 2915 2916 Logically Collective 2917 2918 Input Parameters: 2919 + pre - string to prepend to name or NULL 2920 . oldname - the old, deprecated option 2921 . newname - the new option, or NULL if option is purely removed 2922 . version - a string describing the version of first deprecation, e.g. "3.9" 2923 - info - additional information string, or NULL. 2924 2925 Options Database Keys: 2926 . -options_suppress_deprecated_warnings - do not print deprecation warnings 2927 2928 Notes: 2929 Must be called between PetscOptionsBegin() (or PetscObjectOptionsBegin()) and PetscOptionsEnd(). 2930 Only the proces of rank zero that owns the PetscOptionsItems are argument (managed by PetscOptionsBegin() or 2931 PetscObjectOptionsBegin() prints the information 2932 If newname is provided, the old option is replaced. Otherwise, it remains 2933 in the options database. 2934 If an option is not replaced, the info argument should be used to advise the user 2935 on how to proceed. 2936 There is a limit on the length of the warning printed, so very long strings 2937 provided as info may be truncated. 2938 2939 Level: developer 2940 2941 .seealso: PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsScalar(), PetscOptionsBool(), PetscOptionsString(), PetscOptionsSetValue() 2942 2943 @*/ 2944 PetscErrorCode PetscOptionsDeprecated_Private(PetscOptionItems *PetscOptionsObject,const char oldname[],const char newname[],const char version[],const char info[]) 2945 { 2946 PetscErrorCode ierr; 2947 PetscBool found,quiet; 2948 const char *value; 2949 const char * const quietopt="-options_suppress_deprecated_warnings"; 2950 char msg[4096]; 2951 char *prefix = NULL; 2952 PetscOptions options = NULL; 2953 MPI_Comm comm = PETSC_COMM_SELF; 2954 2955 PetscFunctionBegin; 2956 PetscValidCharPointer(oldname,2); 2957 PetscValidCharPointer(version,4); 2958 if (PetscOptionsObject) { 2959 prefix = PetscOptionsObject->prefix; 2960 options = PetscOptionsObject->options; 2961 comm = PetscOptionsObject->comm; 2962 } 2963 ierr = PetscOptionsFindPair(options,prefix,oldname,&value,&found);CHKERRQ(ierr); 2964 if (found) { 2965 if (newname) { 2966 if (prefix) { 2967 ierr = PetscOptionsPrefixPush(options,prefix);CHKERRQ(ierr); 2968 } 2969 ierr = PetscOptionsSetValue(options,newname,value);CHKERRQ(ierr); 2970 if (prefix) { 2971 ierr = PetscOptionsPrefixPop(options);CHKERRQ(ierr); 2972 } 2973 ierr = PetscOptionsClearValue(options,oldname);CHKERRQ(ierr); 2974 } 2975 quiet = PETSC_FALSE; 2976 ierr = PetscOptionsGetBool(options,NULL,quietopt,&quiet,NULL);CHKERRQ(ierr); 2977 if (!quiet) { 2978 ierr = PetscStrcpy(msg,"** PETSc DEPRECATION WARNING ** : the option ");CHKERRQ(ierr); 2979 ierr = PetscStrcat(msg,oldname);CHKERRQ(ierr); 2980 ierr = PetscStrcat(msg," is deprecated as of version ");CHKERRQ(ierr); 2981 ierr = PetscStrcat(msg,version);CHKERRQ(ierr); 2982 ierr = PetscStrcat(msg," and will be removed in a future release.");CHKERRQ(ierr); 2983 if (newname) { 2984 ierr = PetscStrcat(msg," Please use the option ");CHKERRQ(ierr); 2985 ierr = PetscStrcat(msg,newname);CHKERRQ(ierr); 2986 ierr = PetscStrcat(msg," instead.");CHKERRQ(ierr); 2987 } 2988 if (info) { 2989 ierr = PetscStrcat(msg," ");CHKERRQ(ierr); 2990 ierr = PetscStrcat(msg,info);CHKERRQ(ierr); 2991 } 2992 ierr = PetscStrcat(msg," (Silence this warning with ");CHKERRQ(ierr); 2993 ierr = PetscStrcat(msg,quietopt);CHKERRQ(ierr); 2994 ierr = PetscStrcat(msg,")\n");CHKERRQ(ierr); 2995 ierr = PetscPrintf(comm,msg);CHKERRQ(ierr); 2996 } 2997 } 2998 PetscFunctionReturn(0); 2999 } 3000