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