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