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