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