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