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