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