153acd3b1SBarry Smith /* 23fc1eb6aSBarry Smith Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to 33fc1eb6aSBarry Smith GUI code to display the options and get values from the users. 453acd3b1SBarry Smith 553acd3b1SBarry Smith */ 653acd3b1SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/ 8665c2dedSJed Brown #include <petscviewer.h> 953acd3b1SBarry Smith 102aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None") 112aa6d131SJed Brown 1253acd3b1SBarry Smith /* 1353acd3b1SBarry Smith Keep a linked list of options that have been posted and we are waiting for 143fc1eb6aSBarry Smith user selection. See the manual page for PetscOptionsBegin() 1553acd3b1SBarry Smith 1653acd3b1SBarry Smith Eventually we'll attach this beast to a MPI_Comm 1753acd3b1SBarry Smith */ 18e55864a3SBarry Smith 1953acd3b1SBarry Smith /* 2053acd3b1SBarry Smith Handles setting up the data structure in a call to PetscOptionsBegin() 2153acd3b1SBarry Smith */ 22d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject, MPI_Comm comm, const char prefix[], const char title[], const char mansec[]) 23d71ae5a4SJacob Faibussowitsch { 2453acd3b1SBarry Smith PetscFunctionBegin; 25064a246eSJacob Faibussowitsch if (prefix) PetscValidCharPointer(prefix, 3); 26064a246eSJacob Faibussowitsch PetscValidCharPointer(title, 4); 27064a246eSJacob Faibussowitsch if (mansec) PetscValidCharPointer(mansec, 5); 280eb63584SBarry Smith if (!PetscOptionsObject->alreadyprinted) { 299566063dSJacob Faibussowitsch if (!PetscOptionsHelpPrintedSingleton) PetscCall(PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton, prefix, title, &PetscOptionsObject->alreadyprinted)); 310eb63584SBarry Smith } 3202c9f0b5SLisandro Dalcin PetscOptionsObject->next = NULL; 33e55864a3SBarry Smith PetscOptionsObject->comm = comm; 34e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_FALSE; 35a297a907SKarl Rupp 369566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(prefix, &PetscOptionsObject->prefix)); 379566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(title, &PetscOptionsObject->title)); 3853acd3b1SBarry Smith 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasHelp(PetscOptionsObject->options, &PetscOptionsObject->printhelp)); 40e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) { 4148a46eb9SPierre Jolivet if (!PetscOptionsObject->alreadyprinted) PetscCall((*PetscHelpPrintf)(comm, "----------------------------------------\n%s:\n", title)); 4261b37b28SSatish Balay } 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4453acd3b1SBarry Smith } 4553acd3b1SBarry Smith 463194b578SJed Brown /* 473194b578SJed Brown Handles setting up the data structure in a call to PetscObjectOptionsBegin() 483194b578SJed Brown */ 49d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj, PetscOptionItems *PetscOptionsObject) 50d71ae5a4SJacob Faibussowitsch { 513194b578SJed Brown char title[256]; 523194b578SJed Brown PetscBool flg; 533194b578SJed Brown 543194b578SJed Brown PetscFunctionBegin; 55dbbe0bcdSBarry Smith PetscValidPointer(PetscOptionsObject, 2); 56dbbe0bcdSBarry Smith PetscValidHeader(obj, 1); 57e55864a3SBarry Smith PetscOptionsObject->object = obj; 58e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = obj->optionsprinted; 59a297a907SKarl Rupp 609566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(obj->description, obj->class_name, &flg)); 619566063dSJacob Faibussowitsch if (flg) PetscCall(PetscSNPrintf(title, sizeof(title), "%s options", obj->class_name)); 629566063dSJacob Faibussowitsch else PetscCall(PetscSNPrintf(title, sizeof(title), "%s (%s) options", obj->description, obj->class_name)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsBegin_Private(PetscOptionsObject, obj->comm, obj->prefix, title, obj->mansec)); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 653194b578SJed Brown } 663194b578SJed Brown 6753acd3b1SBarry Smith /* 6853acd3b1SBarry Smith Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd() 6953acd3b1SBarry Smith */ 703ba16761SJacob Faibussowitsch static PetscErrorCode PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscOptionType t, PetscOptionItem *amsopt) 71d71ae5a4SJacob Faibussowitsch { 724416b707SBarry Smith PetscOptionItem next; 733be6e4c3SJed Brown PetscBool valid; 7453acd3b1SBarry Smith 7553acd3b1SBarry Smith PetscFunctionBegin; 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsValidKey(opt, &valid)); 775f80ce2aSJacob Faibussowitsch PetscCheck(valid, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "The option '%s' is not a valid key", opt); 783be6e4c3SJed Brown 799566063dSJacob Faibussowitsch PetscCall(PetscNew(amsopt)); 8002c9f0b5SLisandro Dalcin (*amsopt)->next = NULL; 8153acd3b1SBarry Smith (*amsopt)->set = PETSC_FALSE; 826356e834SBarry Smith (*amsopt)->type = t; 8302c9f0b5SLisandro Dalcin (*amsopt)->data = NULL; 8461b37b28SSatish Balay 859566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(text, &(*amsopt)->text)); 869566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(opt, &(*amsopt)->option)); 879566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(man, &(*amsopt)->man)); 8853acd3b1SBarry Smith 89e55864a3SBarry Smith if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt; 90a297a907SKarl Rupp else { 91e55864a3SBarry Smith next = PetscOptionsObject->next; 9253acd3b1SBarry Smith while (next->next) next = next->next; 9353acd3b1SBarry Smith next->next = *amsopt; 9453acd3b1SBarry Smith } 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9653acd3b1SBarry Smith } 9753acd3b1SBarry Smith 98aee2cecaSBarry Smith /* 993fc1eb6aSBarry Smith PetscScanString - Gets user input via stdin from process and broadcasts to all processes 1003fc1eb6aSBarry Smith 101d083f849SBarry Smith Collective 1023fc1eb6aSBarry Smith 1033fc1eb6aSBarry Smith Input Parameters: 1043fc1eb6aSBarry Smith + commm - communicator for the broadcast, must be PETSC_COMM_WORLD 1053fc1eb6aSBarry Smith . n - length of the string, must be the same on all processes 1063fc1eb6aSBarry Smith - str - location to store input 107aee2cecaSBarry Smith 108aee2cecaSBarry Smith Bugs: 109aee2cecaSBarry Smith . Assumes process 0 of the given communicator has access to stdin 110aee2cecaSBarry Smith 111aee2cecaSBarry Smith */ 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscScanString(MPI_Comm comm, size_t n, char str[]) 113d71ae5a4SJacob Faibussowitsch { 1143fc1eb6aSBarry Smith PetscMPIInt rank, nm; 115aee2cecaSBarry Smith 116aee2cecaSBarry Smith PetscFunctionBegin; 1179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 118dd400576SPatrick Sanan if (rank == 0) { 1195f80ce2aSJacob Faibussowitsch char c = (char)getchar(); 1205f80ce2aSJacob Faibussowitsch size_t i = 0; 1215f80ce2aSJacob Faibussowitsch 122aee2cecaSBarry Smith while (c != '\n' && i < n - 1) { 123aee2cecaSBarry Smith str[i++] = c; 124aee2cecaSBarry Smith c = (char)getchar(); 125aee2cecaSBarry Smith } 126aee2cecaSBarry Smith str[i] = 0; 127aee2cecaSBarry Smith } 1289566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(n, &nm)); 1299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(str, nm, MPI_CHAR, 0, comm)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 131aee2cecaSBarry Smith } 132aee2cecaSBarry Smith 1335b02f95dSBarry Smith /* 1345b02f95dSBarry Smith This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy() 1355b02f95dSBarry Smith */ 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscStrdup(const char s[], char *t[]) 137d71ae5a4SJacob Faibussowitsch { 13802c9f0b5SLisandro Dalcin char *tmp = NULL; 1395b02f95dSBarry Smith 1405b02f95dSBarry Smith PetscFunctionBegin; 1415b02f95dSBarry Smith if (s) { 1425f80ce2aSJacob Faibussowitsch size_t len; 1435f80ce2aSJacob Faibussowitsch 1449566063dSJacob Faibussowitsch PetscCall(PetscStrlen(s, &len)); 1455f80ce2aSJacob Faibussowitsch tmp = (char *)malloc((len + 1) * sizeof(*tmp)); 1465f80ce2aSJacob Faibussowitsch PetscCheck(tmp, PETSC_COMM_SELF, PETSC_ERR_MEM, "No memory to duplicate string"); 147*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(tmp, s, len + 1)); 1485b02f95dSBarry Smith } 1495b02f95dSBarry Smith *t = tmp; 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1515b02f95dSBarry Smith } 1525b02f95dSBarry Smith 153aee2cecaSBarry Smith /* 1543cc1e11dSBarry Smith PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime 155aee2cecaSBarry Smith 15695452b02SPatrick Sanan Notes: 15795452b02SPatrick Sanan this isn't really practical, it is just to demonstrate the principle 158aee2cecaSBarry Smith 1597781c08eSBarry Smith A carriage return indicates no change from the default; but this like -ksp_monitor <stdout> the default is actually not stdout the default 1607781c08eSBarry Smith is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug? 1617781c08eSBarry Smith 162aee2cecaSBarry Smith Bugs: 1637781c08eSBarry Smith + All processes must traverse through the exact same set of option queries due to the call to PetscScanString() 1643cc1e11dSBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 165aee2cecaSBarry Smith - Only works for PetscInt == int, PetscReal == double etc 166aee2cecaSBarry Smith 167811af0c4SBarry Smith Developer Note: 16895452b02SPatrick Sanan Normally the GUI that presents the options the user and retrieves the values would be running in a different 1693cc1e11dSBarry Smith address space and communicating with the PETSc program 1703cc1e11dSBarry Smith 171aee2cecaSBarry Smith */ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject) 173d71ae5a4SJacob Faibussowitsch { 1744416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 1756356e834SBarry Smith char str[512]; 1767781c08eSBarry Smith PetscBool bid; 177a4404d99SBarry Smith PetscReal ir, *valr; 178330cf3c9SBarry Smith PetscInt *vald; 179330cf3c9SBarry Smith size_t i; 1806356e834SBarry Smith 1812a409bb0SBarry Smith PetscFunctionBegin; 1829566063dSJacob Faibussowitsch PetscCall((*PetscPrintf)(PETSC_COMM_WORLD, "%s --------------------\n", PetscOptionsObject->title)); 1836356e834SBarry Smith while (next) { 1846356e834SBarry Smith switch (next->type) { 185d71ae5a4SJacob Faibussowitsch case OPTION_HEAD: 186d71ae5a4SJacob Faibussowitsch break; 187e26ddf31SBarry Smith case OPTION_INT_ARRAY: 1889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1)); 189e26ddf31SBarry Smith vald = (PetscInt *)next->data; 190e26ddf31SBarry Smith for (i = 0; i < next->arraylength; i++) { 1919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT, vald[i])); 19248a46eb9SPierre Jolivet if (i < next->arraylength - 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ",")); 193e26ddf31SBarry Smith } 1949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, ">: %s (%s) ", next->text, next->man)); 1959566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 196e26ddf31SBarry Smith if (str[0]) { 197e26ddf31SBarry Smith PetscToken token; 198e26ddf31SBarry Smith PetscInt n = 0, nmax = next->arraylength, *dvalue = (PetscInt *)next->data, start, end; 199e26ddf31SBarry Smith size_t len; 200e26ddf31SBarry Smith char *value; 201ace3abfcSBarry Smith PetscBool foundrange; 202e26ddf31SBarry Smith 203e26ddf31SBarry Smith next->set = PETSC_TRUE; 204e26ddf31SBarry Smith value = str; 2059566063dSJacob Faibussowitsch PetscCall(PetscTokenCreate(value, ',', &token)); 2069566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token, &value)); 207e26ddf31SBarry Smith while (n < nmax) { 208e26ddf31SBarry Smith if (!value) break; 209e26ddf31SBarry Smith 210e26ddf31SBarry Smith /* look for form d-D where d and D are integers */ 211e26ddf31SBarry Smith foundrange = PETSC_FALSE; 2129566063dSJacob Faibussowitsch PetscCall(PetscStrlen(value, &len)); 213e26ddf31SBarry Smith if (value[0] == '-') i = 2; 214e26ddf31SBarry Smith else i = 1; 215330cf3c9SBarry Smith for (; i < len; i++) { 216e26ddf31SBarry Smith if (value[i] == '-') { 21708401ef6SPierre Jolivet PetscCheck(i != len - 1, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry %s", n, value); 218e26ddf31SBarry Smith value[i] = 0; 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value, &start)); 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value + i + 1, &end)); 22108401ef6SPierre Jolivet PetscCheck(end > start, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry, %s-%s cannot have decreasing list", n, value, value + i + 1); 222cc73adaaSBarry Smith PetscCheck(n + end - start - 1 < nmax, PETSC_COMM_SELF, PETSC_ERR_USER, "Error in %" PetscInt_FMT "-th array entry, not enough space in left in array (%" PetscInt_FMT ") to contain entire range from %" PetscInt_FMT " to %" PetscInt_FMT, n, nmax - n, start, end); 223e26ddf31SBarry Smith for (; start < end; start++) { 2249371c9d4SSatish Balay *dvalue = start; 2259371c9d4SSatish Balay dvalue++; 2269371c9d4SSatish Balay n++; 227e26ddf31SBarry Smith } 228e26ddf31SBarry Smith foundrange = PETSC_TRUE; 229e26ddf31SBarry Smith break; 230e26ddf31SBarry Smith } 231e26ddf31SBarry Smith } 232e26ddf31SBarry Smith if (!foundrange) { 2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value, dvalue)); 234e26ddf31SBarry Smith dvalue++; 235e26ddf31SBarry Smith n++; 236e26ddf31SBarry Smith } 2379566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token, &value)); 238e26ddf31SBarry Smith } 2399566063dSJacob Faibussowitsch PetscCall(PetscTokenDestroy(&token)); 240e26ddf31SBarry Smith } 241e26ddf31SBarry Smith break; 242e26ddf31SBarry Smith case OPTION_REAL_ARRAY: 2439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1)); 244e26ddf31SBarry Smith valr = (PetscReal *)next->data; 245e26ddf31SBarry Smith for (i = 0; i < next->arraylength; i++) { 2469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g", (double)valr[i])); 24748a46eb9SPierre Jolivet if (i < next->arraylength - 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ",")); 248e26ddf31SBarry Smith } 2499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, ">: %s (%s) ", next->text, next->man)); 2509566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 251e26ddf31SBarry Smith if (str[0]) { 252e26ddf31SBarry Smith PetscToken token; 253e26ddf31SBarry Smith PetscInt n = 0, nmax = next->arraylength; 254e26ddf31SBarry Smith PetscReal *dvalue = (PetscReal *)next->data; 255e26ddf31SBarry Smith char *value; 256e26ddf31SBarry Smith 257e26ddf31SBarry Smith next->set = PETSC_TRUE; 258e26ddf31SBarry Smith value = str; 2599566063dSJacob Faibussowitsch PetscCall(PetscTokenCreate(value, ',', &token)); 2609566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token, &value)); 261e26ddf31SBarry Smith while (n < nmax) { 262e26ddf31SBarry Smith if (!value) break; 2639566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToReal(value, dvalue)); 264e26ddf31SBarry Smith dvalue++; 265e26ddf31SBarry Smith n++; 2669566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token, &value)); 267e26ddf31SBarry Smith } 2689566063dSJacob Faibussowitsch PetscCall(PetscTokenDestroy(&token)); 269e26ddf31SBarry Smith } 270e26ddf31SBarry Smith break; 2716356e834SBarry Smith case OPTION_INT: 2729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <%d>: %s (%s) ", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, *(int *)next->data, next->text, next->man)); 2739566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 2743fc1eb6aSBarry Smith if (str[0]) { 275d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 276d25d7f95SJed Brown long long lid; 277d25d7f95SJed Brown sscanf(str, "%lld", &lid); 278cc73adaaSBarry Smith PetscCheck(lid <= PETSC_MAX_INT && lid >= PETSC_MIN_INT, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Argument: -%s%s %lld", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, lid); 279c272547aSJed Brown #else 280d25d7f95SJed Brown long lid; 281d25d7f95SJed Brown sscanf(str, "%ld", &lid); 282cc73adaaSBarry Smith PetscCheck(lid <= PETSC_MAX_INT && lid >= PETSC_MIN_INT, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Argument: -%s%s %ld", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, lid); 283c272547aSJed Brown #endif 284a297a907SKarl Rupp 285d25d7f95SJed Brown next->set = PETSC_TRUE; 286d25d7f95SJed Brown *((PetscInt *)next->data) = (PetscInt)lid; 287aee2cecaSBarry Smith } 288aee2cecaSBarry Smith break; 289aee2cecaSBarry Smith case OPTION_REAL: 2909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <%g>: %s (%s) ", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, *(double *)next->data, next->text, next->man)); 2919566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 2923fc1eb6aSBarry Smith if (str[0]) { 293ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 294a4404d99SBarry Smith sscanf(str, "%e", &ir); 295570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 296570b7f6dSBarry Smith float irtemp; 297570b7f6dSBarry Smith sscanf(str, "%e", &irtemp); 298570b7f6dSBarry Smith ir = irtemp; 299ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE) 300aee2cecaSBarry Smith sscanf(str, "%le", &ir); 301ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 302d9822059SBarry Smith ir = strtoflt128(str, 0); 303d9822059SBarry Smith #else 304513dbe71SLisandro Dalcin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Unknown scalar type"); 305a4404d99SBarry Smith #endif 306aee2cecaSBarry Smith next->set = PETSC_TRUE; 307aee2cecaSBarry Smith *((PetscReal *)next->data) = ir; 308aee2cecaSBarry Smith } 309aee2cecaSBarry Smith break; 3107781c08eSBarry Smith case OPTION_BOOL: 3119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <%s>: %s (%s) ", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, *(PetscBool *)next->data ? "true" : "false", next->text, next->man)); 3129566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 3137781c08eSBarry Smith if (str[0]) { 3149566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToBool(str, &bid)); 3157781c08eSBarry Smith next->set = PETSC_TRUE; 3167781c08eSBarry Smith *((PetscBool *)next->data) = bid; 3177781c08eSBarry Smith } 3187781c08eSBarry Smith break; 319aee2cecaSBarry Smith case OPTION_STRING: 3209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "-%s%s <%s>: %s (%s) ", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", next->option + 1, (char *)next->data, next->text, next->man)); 3219566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 3223fc1eb6aSBarry Smith if (str[0]) { 323aee2cecaSBarry Smith next->set = PETSC_TRUE; 32464facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3259566063dSJacob Faibussowitsch PetscCall(PetscStrdup(str, (char **)&next->data)); 3266356e834SBarry Smith } 3276356e834SBarry Smith break; 328a264d7a6SBarry Smith case OPTION_FLIST: 3299566063dSJacob Faibussowitsch PetscCall(PetscFunctionListPrintTypes(PETSC_COMM_WORLD, stdout, PetscOptionsObject->prefix, next->option, next->text, next->man, next->flist, (char *)next->data, (char *)next->data)); 3309566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD, 512, str)); 3313cc1e11dSBarry Smith if (str[0]) { 332e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_TRUE; 3333cc1e11dSBarry Smith next->set = PETSC_TRUE; 33464facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3359566063dSJacob Faibussowitsch PetscCall(PetscStrdup(str, (char **)&next->data)); 3363cc1e11dSBarry Smith } 3373cc1e11dSBarry Smith break; 338d71ae5a4SJacob Faibussowitsch default: 339d71ae5a4SJacob Faibussowitsch break; 3406356e834SBarry Smith } 3416356e834SBarry Smith next = next->next; 3426356e834SBarry Smith } 3433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3446356e834SBarry Smith } 3456356e834SBarry Smith 346e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 347e04113cfSBarry Smith #include <petscviewersaws.h> 348d5649816SBarry Smith 349d5649816SBarry Smith static int count = 0; 350d5649816SBarry Smith 351d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSAWsDestroy(void) 352d71ae5a4SJacob Faibussowitsch { 3532657e9d9SBarry Smith PetscFunctionBegin; 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355d5649816SBarry Smith } 356d5649816SBarry Smith 3579c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n" 358a8d69d7bSBarry Smith "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n" 359a8d69d7bSBarry Smith "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n" 360d1fc0251SBarry Smith "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n" 36164bbc9efSBarry Smith "<script>\n" 36264bbc9efSBarry Smith "jQuery(document).ready(function() {\n" 36364bbc9efSBarry Smith "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n" 36464bbc9efSBarry Smith "})\n" 36564bbc9efSBarry Smith "</script>\n" 36664bbc9efSBarry Smith "</head>\n"; 3671423471aSBarry Smith 3681423471aSBarry Smith /* Determines the size and style of the scroll region where PETSc options selectable from users are displayed */ 3691423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>"; 37064bbc9efSBarry Smith 371b3506946SBarry Smith /* 3727781c08eSBarry Smith PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs 373b3506946SBarry Smith 374b3506946SBarry Smith Bugs: 375b3506946SBarry Smith + All processes must traverse through the exact same set of option queries do to the call to PetscScanString() 376b3506946SBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 377b3506946SBarry Smith - Only works for PetscInt == int, PetscReal == double etc 378b3506946SBarry Smith 379b3506946SBarry Smith */ 380d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject) 381d71ae5a4SJacob Faibussowitsch { 3824416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 383d5649816SBarry Smith static int mancount = 0; 384b3506946SBarry Smith char options[16]; 3857aab2a10SBarry Smith PetscBool changedmethod = PETSC_FALSE; 386a23eb982SSurtai Han PetscBool stopasking = PETSC_FALSE; 38788a9752cSBarry Smith char manname[16], textname[16]; 3882657e9d9SBarry Smith char dir[1024]; 389b3506946SBarry Smith 3902a409bb0SBarry Smith PetscFunctionBegin; 391b3506946SBarry Smith /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */ 392a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(options, PETSC_STATIC_ARRAY_LENGTH(options), "Options_%d", count++)); 393a297a907SKarl Rupp 3947325c4a4SBarry Smith PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */ 3951bc75a8dSBarry Smith 3969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", "_title")); 397792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &PetscOptionsObject->title, 1, SAWs_READ, SAWs_STRING)); 3989566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", "prefix")); 399792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &PetscOptionsObject->pprefix, 1, SAWs_READ, SAWs_STRING)); 400792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, ("/PETSc/Options/ChangedMethod", &changedmethod, 1, SAWs_WRITE, SAWs_BOOLEAN)); 401792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, ("/PETSc/Options/StopAsking", &stopasking, 1, SAWs_WRITE, SAWs_BOOLEAN)); 402b3506946SBarry Smith 403b3506946SBarry Smith while (next) { 404a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(manname, PETSC_STATIC_ARRAY_LENGTH(manname), "_man_%d", mancount)); 4059566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", manname)); 406792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &next->man, 1, SAWs_READ, SAWs_STRING)); 407a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(textname, PETSC_STATIC_ARRAY_LENGTH(textname), "_text_%d", mancount++)); 4089566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", textname)); 409792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &next->text, 1, SAWs_READ, SAWs_STRING)); 4109f32e415SBarry Smith 411b3506946SBarry Smith switch (next->type) { 412d71ae5a4SJacob Faibussowitsch case OPTION_HEAD: 413d71ae5a4SJacob Faibussowitsch break; 414b3506946SBarry Smith case OPTION_INT_ARRAY: 4159566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 416792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, next->arraylength, SAWs_WRITE, SAWs_INT)); 417b3506946SBarry Smith break; 418b3506946SBarry Smith case OPTION_REAL_ARRAY: 4199566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 420792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, next->arraylength, SAWs_WRITE, SAWs_DOUBLE)); 421b3506946SBarry Smith break; 422b3506946SBarry Smith case OPTION_INT: 4239566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 424792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, 1, SAWs_WRITE, SAWs_INT)); 425b3506946SBarry Smith break; 426b3506946SBarry Smith case OPTION_REAL: 4279566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 428792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, 1, SAWs_WRITE, SAWs_DOUBLE)); 429b3506946SBarry Smith break; 4307781c08eSBarry Smith case OPTION_BOOL: 4319566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 432792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, 1, SAWs_WRITE, SAWs_BOOLEAN)); 4331ae3d29cSBarry Smith break; 4347781c08eSBarry Smith case OPTION_BOOL_ARRAY: 4359566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 436792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, next->arraylength, SAWs_WRITE, SAWs_BOOLEAN)); 43771f08665SBarry Smith break; 438b3506946SBarry Smith case OPTION_STRING: 4399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 440792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &next->data, 1, SAWs_WRITE, SAWs_STRING)); 4411ae3d29cSBarry Smith break; 4421ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 4439566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 444792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, next->data, next->arraylength, SAWs_WRITE, SAWs_STRING)); 445b3506946SBarry Smith break; 4469371c9d4SSatish Balay case OPTION_FLIST: { 447a264d7a6SBarry Smith PetscInt ntext; 4489566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 449792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &next->data, 1, SAWs_WRITE, SAWs_STRING)); 4509566063dSJacob Faibussowitsch PetscCall(PetscFunctionListGet(next->flist, (const char ***)&next->edata, &ntext)); 451792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Legal_Variable_Values, (dir, ntext, next->edata)); 4529371c9d4SSatish Balay } break; 4539371c9d4SSatish Balay case OPTION_ELIST: { 454a264d7a6SBarry Smith PetscInt ntext = next->nlist; 4559566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 456792fecdfSBarry Smith PetscCallSAWs(SAWs_Register, (dir, &next->data, 1, SAWs_WRITE, SAWs_STRING)); 4579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((ntext + 1), (char ***)&next->edata)); 4589566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(next->edata, next->list, ntext * sizeof(char *))); 459792fecdfSBarry Smith PetscCallSAWs(SAWs_Set_Legal_Variable_Values, (dir, ntext, next->edata)); 4609371c9d4SSatish Balay } break; 461d71ae5a4SJacob Faibussowitsch default: 462d71ae5a4SJacob Faibussowitsch break; 463b3506946SBarry Smith } 464b3506946SBarry Smith next = next->next; 465b3506946SBarry Smith } 466b3506946SBarry Smith 467b3506946SBarry Smith /* wait until accessor has unlocked the memory */ 468792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Header, ("index.html", OptionsHeader)); 469792fecdfSBarry Smith PetscCallSAWs(SAWs_Push_Body, ("index.html", 2, OptionsBodyBottom)); 4709566063dSJacob Faibussowitsch PetscCall(PetscSAWsBlock()); 471792fecdfSBarry Smith PetscCallSAWs(SAWs_Pop_Header, ("index.html")); 472792fecdfSBarry Smith PetscCallSAWs(SAWs_Pop_Body, ("index.html", 2)); 473b3506946SBarry Smith 47488a9752cSBarry Smith /* determine if any values have been set in GUI */ 47583355fc5SBarry Smith next = PetscOptionsObject->next; 47688a9752cSBarry Smith while (next) { 4779566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Options/%s", next->option)); 478792fecdfSBarry Smith PetscCallSAWs(SAWs_Selected, (dir, (int *)&next->set)); 47988a9752cSBarry Smith next = next->next; 48088a9752cSBarry Smith } 48188a9752cSBarry Smith 482b3506946SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 483f7b25cf6SBarry Smith if (changedmethod) PetscOptionsObject->count = -2; 484b3506946SBarry Smith 485a23eb982SSurtai Han if (stopasking) { 486a23eb982SSurtai Han PetscOptionsPublish = PETSC_FALSE; 48712655325SBarry Smith PetscOptionsObject->count = 0; //do not ask for same thing again 488a23eb982SSurtai Han } 489a23eb982SSurtai Han 490792fecdfSBarry Smith PetscCallSAWs(SAWs_Delete, ("/PETSc/Options")); 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 492b3506946SBarry Smith } 493b3506946SBarry Smith #endif 494b3506946SBarry Smith 495d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject) 496d71ae5a4SJacob Faibussowitsch { 4974416b707SBarry Smith PetscOptionItem last; 4986356e834SBarry Smith char option[256], value[1024], tmp[32]; 499330cf3c9SBarry Smith size_t j; 50053acd3b1SBarry Smith 50153acd3b1SBarry Smith PetscFunctionBegin; 50283355fc5SBarry Smith if (PetscOptionsObject->next) { 50383355fc5SBarry Smith if (!PetscOptionsObject->count) { 504a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS) 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsSAWsInput(PetscOptionsObject)); 506b3506946SBarry Smith #else 5079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetFromTextInput(PetscOptionsObject)); 508b3506946SBarry Smith #endif 509aee2cecaSBarry Smith } 510aee2cecaSBarry Smith } 5116356e834SBarry Smith 5129566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->title)); 5136356e834SBarry Smith 514e26ddf31SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 515e55864a3SBarry Smith if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2; 5167a72a596SBarry Smith /* reset alreadyprinted flag */ 517e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = PETSC_FALSE; 518e55864a3SBarry Smith if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE; 519e55864a3SBarry Smith PetscOptionsObject->object = NULL; 52053acd3b1SBarry Smith 521e55864a3SBarry Smith while (PetscOptionsObject->next) { 522e55864a3SBarry Smith if (PetscOptionsObject->next->set) { 523e55864a3SBarry Smith if (PetscOptionsObject->prefix) { 524*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(option, "-", sizeof(option))); 525*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(option, PetscOptionsObject->prefix, sizeof(option))); 526*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(option, PetscOptionsObject->next->option + 1, sizeof(option))); 527*c6a7a370SJeremy L Thompson } else PetscCall(PetscStrncpy(option, PetscOptionsObject->next->option, sizeof(option))); 5286356e834SBarry Smith 529e55864a3SBarry Smith switch (PetscOptionsObject->next->type) { 530d71ae5a4SJacob Faibussowitsch case OPTION_HEAD: 531d71ae5a4SJacob Faibussowitsch break; 5326356e834SBarry Smith case OPTION_INT_ARRAY: 533a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%d", (int)((PetscInt *)PetscOptionsObject->next->data)[0])); 534e55864a3SBarry Smith for (j = 1; j < PetscOptionsObject->next->arraylength; j++) { 535a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(tmp, PETSC_STATIC_ARRAY_LENGTH(tmp), "%d", (int)((PetscInt *)PetscOptionsObject->next->data)[j])); 536*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, ",", sizeof(value))); 537*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, tmp, sizeof(value))); 5386356e834SBarry Smith } 5396356e834SBarry Smith break; 540d71ae5a4SJacob Faibussowitsch case OPTION_INT: 541a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%d", (int)*(PetscInt *)PetscOptionsObject->next->data)); 542d71ae5a4SJacob Faibussowitsch break; 543d71ae5a4SJacob Faibussowitsch case OPTION_REAL: 544a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%g", (double)*(PetscReal *)PetscOptionsObject->next->data)); 545d71ae5a4SJacob Faibussowitsch break; 5466356e834SBarry Smith case OPTION_REAL_ARRAY: 547a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%g", (double)((PetscReal *)PetscOptionsObject->next->data)[0])); 548e55864a3SBarry Smith for (j = 1; j < PetscOptionsObject->next->arraylength; j++) { 549a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(tmp, PETSC_STATIC_ARRAY_LENGTH(tmp), "%g", (double)((PetscReal *)PetscOptionsObject->next->data)[j])); 550*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, ",", sizeof(value))); 551*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, tmp, sizeof(value))); 5526356e834SBarry Smith } 5536356e834SBarry Smith break; 554050cccc3SHong Zhang case OPTION_SCALAR_ARRAY: 555a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%g+%gi", (double)PetscRealPart(((PetscScalar *)PetscOptionsObject->next->data)[0]), (double)PetscImaginaryPart(((PetscScalar *)PetscOptionsObject->next->data)[0]))); 556050cccc3SHong Zhang for (j = 1; j < PetscOptionsObject->next->arraylength; j++) { 557a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(tmp, PETSC_STATIC_ARRAY_LENGTH(tmp), "%g+%gi", (double)PetscRealPart(((PetscScalar *)PetscOptionsObject->next->data)[j]), (double)PetscImaginaryPart(((PetscScalar *)PetscOptionsObject->next->data)[j]))); 558*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, ",", sizeof(value))); 559*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, tmp, sizeof(value))); 560050cccc3SHong Zhang } 561050cccc3SHong Zhang break; 562d71ae5a4SJacob Faibussowitsch case OPTION_BOOL: 563a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%d", *(int *)PetscOptionsObject->next->data)); 564d71ae5a4SJacob Faibussowitsch break; 5657781c08eSBarry Smith case OPTION_BOOL_ARRAY: 566a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%d", (int)((PetscBool *)PetscOptionsObject->next->data)[0])); 567e55864a3SBarry Smith for (j = 1; j < PetscOptionsObject->next->arraylength; j++) { 568a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(tmp, PETSC_STATIC_ARRAY_LENGTH(tmp), "%d", (int)((PetscBool *)PetscOptionsObject->next->data)[j])); 569*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, ",", sizeof(value))); 570*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, tmp, sizeof(value))); 5711ae3d29cSBarry Smith } 5721ae3d29cSBarry Smith break; 573d71ae5a4SJacob Faibussowitsch case OPTION_FLIST: 574*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(value, (char *)PetscOptionsObject->next->data, sizeof(value))); 575d71ae5a4SJacob Faibussowitsch break; 576d71ae5a4SJacob Faibussowitsch case OPTION_ELIST: 577*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(value, (char *)PetscOptionsObject->next->data, sizeof(value))); 578d71ae5a4SJacob Faibussowitsch break; 579d71ae5a4SJacob Faibussowitsch case OPTION_STRING: 580*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(value, (char *)PetscOptionsObject->next->data, sizeof(value))); 581d71ae5a4SJacob Faibussowitsch break; 5821ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 583a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, PETSC_STATIC_ARRAY_LENGTH(value), "%s", ((char **)PetscOptionsObject->next->data)[0])); 584e55864a3SBarry Smith for (j = 1; j < PetscOptionsObject->next->arraylength; j++) { 585a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(tmp, PETSC_STATIC_ARRAY_LENGTH(tmp), "%s", ((char **)PetscOptionsObject->next->data)[j])); 586*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, ",", sizeof(value))); 587*c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(value, tmp, sizeof(value))); 5881ae3d29cSBarry Smith } 5896356e834SBarry Smith break; 5906356e834SBarry Smith } 5919566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(PetscOptionsObject->options, option, value)); 5926356e834SBarry Smith } 59348a46eb9SPierre Jolivet if (PetscOptionsObject->next->type == OPTION_ELIST) PetscCall(PetscStrNArrayDestroy(PetscOptionsObject->next->nlist, (char ***)&PetscOptionsObject->next->list)); 5949566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->text)); 5959566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->option)); 5969566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->man)); 5979566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->edata)); 598c979a496SBarry Smith 59983355fc5SBarry Smith if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)) { 60083355fc5SBarry Smith free(PetscOptionsObject->next->data); 601c979a496SBarry Smith } else { 6029566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->data)); 603c979a496SBarry Smith } 6047781c08eSBarry Smith 60583355fc5SBarry Smith last = PetscOptionsObject->next; 60683355fc5SBarry Smith PetscOptionsObject->next = PetscOptionsObject->next->next; 6079566063dSJacob Faibussowitsch PetscCall(PetscFree(last)); 6086356e834SBarry Smith } 6099566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->prefix)); 61002c9f0b5SLisandro Dalcin PetscOptionsObject->next = NULL; 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61253acd3b1SBarry Smith } 61353acd3b1SBarry Smith 61488aa4217SBarry Smith /*MC 61553acd3b1SBarry Smith PetscOptionsEnum - Gets the enum value for a particular option in the database. 61653acd3b1SBarry Smith 617811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 61853acd3b1SBarry Smith 61988aa4217SBarry Smith Synopsis: 62088aa4217SBarry Smith #include "petscsys.h" 6213a89f35bSSatish Balay PetscErrorCode PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool *set) 62288aa4217SBarry Smith 62353acd3b1SBarry Smith Input Parameters: 62453acd3b1SBarry Smith + opt - option name 62553acd3b1SBarry Smith . text - short string that describes the option 62653acd3b1SBarry Smith . man - manual page with additional information on option 62753acd3b1SBarry Smith . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 6280fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6290fdccdaeSBarry Smith $ PetscOptionsEnum(..., obj->value,&object->value,...) or 6300fdccdaeSBarry Smith $ value = defaultvalue 6310fdccdaeSBarry Smith $ PetscOptionsEnum(..., value,&value,&flg); 6320fdccdaeSBarry Smith $ if (flg) { 63353acd3b1SBarry Smith 634d8d19677SJose E. Roman Output Parameters: 63553acd3b1SBarry Smith + value - the value to return 636811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 63753acd3b1SBarry Smith 63853acd3b1SBarry Smith Level: beginner 63953acd3b1SBarry Smith 64095452b02SPatrick Sanan Notes: 641811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 64253acd3b1SBarry Smith 643811af0c4SBarry Smith list is usually something like `PCASMTypes` or some other predefined list of enum names 64453acd3b1SBarry Smith 6452efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 6462efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 6472efd9cb1SBarry Smith 648989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 649989712b9SBarry Smith 650db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 651db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 652db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 653db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 654c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 655db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 656db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 65788aa4217SBarry Smith M*/ 65888aa4217SBarry Smith 659d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], const char *const *list, PetscEnum currentvalue, PetscEnum *value, PetscBool *set) 660d71ae5a4SJacob Faibussowitsch { 66153acd3b1SBarry Smith PetscInt ntext = 0; 662aa5bb8c0SSatish Balay PetscInt tval; 663ace3abfcSBarry Smith PetscBool tflg; 66453acd3b1SBarry Smith 66553acd3b1SBarry Smith PetscFunctionBegin; 666ad540459SPierre Jolivet while (list[ntext++]) PetscCheck(ntext <= 50, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument appears to be wrong or have more than 50 entries"); 66708401ef6SPierre Jolivet PetscCheck(ntext >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix"); 66853acd3b1SBarry Smith ntext -= 3; 6699566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList_Private(PetscOptionsObject, opt, text, man, list, ntext, list[currentvalue], &tval, &tflg)); 670aa5bb8c0SSatish Balay /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 671aa5bb8c0SSatish Balay if (tflg) *value = (PetscEnum)tval; 672aa5bb8c0SSatish Balay if (set) *set = tflg; 6733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67453acd3b1SBarry Smith } 67553acd3b1SBarry Smith 67688aa4217SBarry Smith /*MC 677d3e47460SLisandro Dalcin PetscOptionsEnumArray - Gets an array of enum values for a particular 678d3e47460SLisandro Dalcin option in the database. 679d3e47460SLisandro Dalcin 680811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 681d3e47460SLisandro Dalcin 68288aa4217SBarry Smith Synopsis: 68388aa4217SBarry Smith #include "petscsys.h" 6843a89f35bSSatish Balay PetscErrorCode PetscOptionsEnumArray(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool *set) 68588aa4217SBarry Smith 686d3e47460SLisandro Dalcin Input Parameters: 687d3e47460SLisandro Dalcin + opt - the option one is seeking 688d3e47460SLisandro Dalcin . text - short string describing option 689d3e47460SLisandro Dalcin . man - manual page for option 69022d58ec6SMatthew G. Knepley . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 691c89788bdSBarry Smith - n - maximum number of values allowed in the value array 692d3e47460SLisandro Dalcin 693d8d19677SJose E. Roman Output Parameters: 694d3e47460SLisandro Dalcin + value - location to copy values 695d3e47460SLisandro Dalcin . n - actual number of values found 696811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 697d3e47460SLisandro Dalcin 698d3e47460SLisandro Dalcin Level: beginner 699d3e47460SLisandro Dalcin 700d3e47460SLisandro Dalcin Notes: 701d3e47460SLisandro Dalcin The array must be passed as a comma separated list. 702d3e47460SLisandro Dalcin 703d3e47460SLisandro Dalcin There must be no intervening spaces between the values. 704d3e47460SLisandro Dalcin 705811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 706d3e47460SLisandro Dalcin 707db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 708db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 709db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 710c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 711db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 712db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRealArray()` 71388aa4217SBarry Smith M*/ 71488aa4217SBarry Smith 715d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], const char *const *list, PetscEnum value[], PetscInt *n, PetscBool *set) 716d71ae5a4SJacob Faibussowitsch { 717d3e47460SLisandro Dalcin PetscInt i, nlist = 0; 7184416b707SBarry Smith PetscOptionItem amsopt; 719d3e47460SLisandro Dalcin 720d3e47460SLisandro Dalcin PetscFunctionBegin; 72108401ef6SPierre Jolivet while (list[nlist++]) PetscCheck(nlist <= 50, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument appears to be wrong or have more than 50 entries"); 72208401ef6SPierre Jolivet PetscCheck(nlist >= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "List argument must have at least two entries: typename and type prefix"); 723d3e47460SLisandro Dalcin nlist -= 3; /* drop enum name, prefix, and null termination */ 724d3e47460SLisandro Dalcin if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */ 725d3e47460SLisandro Dalcin PetscEnum *vals; 7269566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_INT_ARRAY /*XXX OPTION_ENUM_ARRAY*/, &amsopt)); 7279566063dSJacob Faibussowitsch PetscCall(PetscStrNArrayallocpy(nlist, list, (char ***)&amsopt->list)); 728d3e47460SLisandro Dalcin amsopt->nlist = nlist; 7299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n, (PetscEnum **)&amsopt->data)); 730d3e47460SLisandro Dalcin amsopt->arraylength = *n; 731d3e47460SLisandro Dalcin vals = (PetscEnum *)amsopt->data; 732d3e47460SLisandro Dalcin for (i = 0; i < *n; i++) vals[i] = value[i]; 733d3e47460SLisandro Dalcin } 7349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, list, value, n, set)); 735d3e47460SLisandro Dalcin if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 7369566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%s", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, list[value[0]])); 7379566063dSJacob Faibussowitsch for (i = 1; i < *n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ",%s", list[value[i]])); 7389566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ">: %s (choose from)", text)); 7399566063dSJacob Faibussowitsch for (i = 0; i < nlist; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " %s", list[i])); 7409566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " (%s)\n", ManSection(man))); 741d3e47460SLisandro Dalcin } 7423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 743d3e47460SLisandro Dalcin } 744d3e47460SLisandro Dalcin 74588aa4217SBarry Smith /*MC 7465a856986SBarry Smith PetscOptionsBoundedInt - Gets an integer value greater than or equal a given bound for a particular option in the database. 7475a856986SBarry Smith 748811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 7495a856986SBarry Smith 75088aa4217SBarry Smith Synopsis: 75188aa4217SBarry Smith #include "petscsys.h" 75269d47153SPierre Jolivet PetscErrorCode PetscOptionsBoundedInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg,PetscInt bound) 75388aa4217SBarry Smith 7545a856986SBarry Smith Input Parameters: 7555a856986SBarry Smith + opt - option name 7565a856986SBarry Smith . text - short string that describes the option 7575a856986SBarry Smith . man - manual page with additional information on option 7585a856986SBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 75969d47153SPierre Jolivet .vb 76069d47153SPierre Jolivet PetscOptionsInt(..., obj->value,&obj->value,...) 76169d47153SPierre Jolivet .ve 76269d47153SPierre Jolivet or 76369d47153SPierre Jolivet .vb 76469d47153SPierre Jolivet value = defaultvalue 76569d47153SPierre Jolivet PetscOptionsInt(..., value,&value,&flg); 76669d47153SPierre Jolivet if (flg) { 76769d47153SPierre Jolivet .ve 76888aa4217SBarry Smith - bound - the requested value should be greater than or equal this bound or an error is generated 7695a856986SBarry Smith 770d8d19677SJose E. Roman Output Parameters: 7715a856986SBarry Smith + value - the integer value to return 772811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 7735a856986SBarry Smith 7745a856986SBarry Smith Notes: 7755a856986SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 7765a856986SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 7775a856986SBarry Smith 7785a856986SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 7795a856986SBarry Smith 780811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 7815a856986SBarry Smith 7825a856986SBarry Smith Level: beginner 7835a856986SBarry Smith 784db781477SPatrick Sanan .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 785db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 786db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 787db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 788c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 789db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 790db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 79188aa4217SBarry Smith M*/ 7925a856986SBarry Smith 79388aa4217SBarry Smith /*MC 7945a856986SBarry Smith PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database. 7955a856986SBarry Smith 796811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 7975a856986SBarry Smith 79888aa4217SBarry Smith Synopsis: 79988aa4217SBarry Smith #include "petscsys.h" 8003a89f35bSSatish Balay PetscErrorCode PetscOptionsRangeInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg,PetscInt lb,PetscInt ub) 80188aa4217SBarry Smith 8025a856986SBarry Smith Input Parameters: 8035a856986SBarry Smith + opt - option name 8045a856986SBarry Smith . text - short string that describes the option 8055a856986SBarry Smith . man - manual page with additional information on option 8065a856986SBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 80788aa4217SBarry Smith $ PetscOptionsInt(..., obj->value,&obj->value,...) or 8085a856986SBarry Smith $ value = defaultvalue 8095a856986SBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 8105a856986SBarry Smith $ if (flg) { 81188aa4217SBarry Smith . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 81288aa4217SBarry Smith - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 8135a856986SBarry Smith 814d8d19677SJose E. Roman Output Parameters: 8155a856986SBarry Smith + value - the integer value to return 816811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 8175a856986SBarry Smith 8185a856986SBarry Smith Notes: 8195a856986SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 8205a856986SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 8215a856986SBarry Smith 8225a856986SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 8235a856986SBarry Smith 824811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 8255a856986SBarry Smith 8265a856986SBarry Smith Level: beginner 8275a856986SBarry Smith 828db781477SPatrick Sanan .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 829db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 830db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 831db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 832c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 833db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 834db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 83588aa4217SBarry Smith M*/ 8365a856986SBarry Smith 83788aa4217SBarry Smith /*MC 83853acd3b1SBarry Smith PetscOptionsInt - Gets the integer value for a particular option in the database. 83953acd3b1SBarry Smith 840811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 84153acd3b1SBarry Smith 84288aa4217SBarry Smith Synopsis: 84388aa4217SBarry Smith #include "petscsys.h" 8446eeb591dSVaclav Hapla PetscErrorCode PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg)) 84588aa4217SBarry Smith 84653acd3b1SBarry Smith Input Parameters: 84753acd3b1SBarry Smith + opt - option name 84853acd3b1SBarry Smith . text - short string that describes the option 84953acd3b1SBarry Smith . man - manual page with additional information on option 8500fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 85188aa4217SBarry Smith $ PetscOptionsInt(..., obj->value,&obj->value,...) or 8520fdccdaeSBarry Smith $ value = defaultvalue 8530fdccdaeSBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 8540fdccdaeSBarry Smith $ if (flg) { 85553acd3b1SBarry Smith 856d8d19677SJose E. Roman Output Parameters: 85753acd3b1SBarry Smith + value - the integer value to return 858811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 85953acd3b1SBarry Smith 86095452b02SPatrick Sanan Notes: 86195452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 8622efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 8632efd9cb1SBarry Smith 864989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 865989712b9SBarry Smith 866811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 86753acd3b1SBarry Smith 8685a856986SBarry Smith Level: beginner 8695a856986SBarry Smith 870db781477SPatrick Sanan .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 871db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 872db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 873db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 874c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 875db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 876db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 87788aa4217SBarry Smith M*/ 87888aa4217SBarry Smith 879d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscInt currentvalue, PetscInt *value, PetscBool *set, PetscInt lb, PetscInt ub) 880d71ae5a4SJacob Faibussowitsch { 8814416b707SBarry Smith PetscOptionItem amsopt; 88212655325SBarry Smith PetscBool wasset; 88353acd3b1SBarry Smith 88453acd3b1SBarry Smith PetscFunctionBegin; 88508401ef6SPierre Jolivet PetscCheck(currentvalue >= lb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Current value %" PetscInt_FMT " less than allowed bound %" PetscInt_FMT, currentvalue, lb); 88608401ef6SPierre Jolivet PetscCheck(currentvalue <= ub, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Current value %" PetscInt_FMT " greater than allowed bound %" PetscInt_FMT, currentvalue, ub); 887e55864a3SBarry Smith if (!PetscOptionsObject->count) { 8889566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_INT, &amsopt)); 8899566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscInt), &amsopt->data)); 89012655325SBarry Smith *(PetscInt *)amsopt->data = currentvalue; 8913e211508SBarry Smith 8929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, ¤tvalue, &wasset)); 893ad540459SPierre Jolivet if (wasset) *(PetscInt *)amsopt->data = currentvalue; 894af6d86caSBarry Smith } 8959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, &wasset)); 896cc73adaaSBarry Smith PetscCheck(!wasset || *value >= lb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Newly set value %" PetscInt_FMT " less than allowed bound %" PetscInt_FMT, *value, lb); 897cc73adaaSBarry Smith PetscCheck(!wasset || *value <= ub, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Newly set value %" PetscInt_FMT " greater than allowed bound %" PetscInt_FMT, *value, ub); 89844ef3d73SBarry Smith if (set) *set = wasset; 899e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9009566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <now %" PetscInt_FMT " : formerly %" PetscInt_FMT ">: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, wasset && value ? *value : currentvalue, currentvalue, text, ManSection(man))); 90153acd3b1SBarry Smith } 9023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90353acd3b1SBarry Smith } 90453acd3b1SBarry Smith 90588aa4217SBarry Smith /*MC 90653acd3b1SBarry Smith PetscOptionsString - Gets the string value for a particular option in the database. 90753acd3b1SBarry Smith 908811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 90953acd3b1SBarry Smith 91088aa4217SBarry Smith Synopsis: 91188aa4217SBarry Smith #include "petscsys.h" 9123a89f35bSSatish Balay PetscErrorCode PetscOptionsString(const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool *set) 91388aa4217SBarry Smith 91453acd3b1SBarry Smith Input Parameters: 91553acd3b1SBarry Smith + opt - option name 91653acd3b1SBarry Smith . text - short string that describes the option 91753acd3b1SBarry Smith . man - manual page with additional information on option 9180fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 919bcbf2dc5SJed Brown - len - length of the result string including null terminator 92053acd3b1SBarry Smith 921d8d19677SJose E. Roman Output Parameters: 92253acd3b1SBarry Smith + value - the value to return 923811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 92453acd3b1SBarry Smith 92553acd3b1SBarry Smith Level: beginner 92653acd3b1SBarry Smith 92795452b02SPatrick Sanan Notes: 928811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 92953acd3b1SBarry Smith 9307fccdfe4SBarry Smith 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). 9317fccdfe4SBarry Smith 9322efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 9332efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 9342efd9cb1SBarry Smith 935989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 936989712b9SBarry Smith 937db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 938db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 939db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 940db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 941c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 942db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 943db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 94488aa4217SBarry Smith M*/ 94588aa4217SBarry Smith 946d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsString_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], const char currentvalue[], char value[], size_t len, PetscBool *set) 947d71ae5a4SJacob Faibussowitsch { 9484416b707SBarry Smith PetscOptionItem amsopt; 94944ef3d73SBarry Smith PetscBool lset; 95053acd3b1SBarry Smith 95153acd3b1SBarry Smith PetscFunctionBegin; 9521a1499c8SBarry Smith if (!PetscOptionsObject->count) { 9539566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_STRING, &amsopt)); 95464facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 9559566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "", (char **)&amsopt->data)); 956af6d86caSBarry Smith } 9579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, len, &lset)); 95844ef3d73SBarry Smith if (set) *set = lset; 959e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9609566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <now %s : formerly %s>: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, lset && value ? value : currentvalue, currentvalue, text, ManSection(man))); 96153acd3b1SBarry Smith } 9623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96353acd3b1SBarry Smith } 96453acd3b1SBarry Smith 96588aa4217SBarry Smith /*MC 966811af0c4SBarry Smith PetscOptionsReal - Gets the `PetscReal` value for a particular option in the database. 96753acd3b1SBarry Smith 968811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 96953acd3b1SBarry Smith 97088aa4217SBarry Smith Synopsis: 97188aa4217SBarry Smith #include "petscsys.h" 9723a89f35bSSatish Balay PetscErrorCode PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 97388aa4217SBarry Smith 97453acd3b1SBarry Smith Input Parameters: 97553acd3b1SBarry Smith + opt - option name 97653acd3b1SBarry Smith . text - short string that describes the option 97753acd3b1SBarry Smith . man - manual page with additional information on option 9780fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 97988aa4217SBarry Smith $ PetscOptionsReal(..., obj->value,&obj->value,...) or 9800fdccdaeSBarry Smith $ value = defaultvalue 9810fdccdaeSBarry Smith $ PetscOptionsReal(..., value,&value,&flg); 9820fdccdaeSBarry Smith $ if (flg) { 98353acd3b1SBarry Smith 984d8d19677SJose E. Roman Output Parameters: 98553acd3b1SBarry Smith + value - the value to return 986811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 98753acd3b1SBarry Smith 98895452b02SPatrick Sanan Notes: 98995452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 9902efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 9912efd9cb1SBarry Smith 992989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 993989712b9SBarry Smith 994811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 99553acd3b1SBarry Smith 9965a856986SBarry Smith Level: beginner 9975a856986SBarry Smith 998db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 999db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1000db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1001db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1002c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1003db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1004db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 100588aa4217SBarry Smith M*/ 100688aa4217SBarry Smith 1007d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscReal currentvalue, PetscReal *value, PetscBool *set) 1008d71ae5a4SJacob Faibussowitsch { 10094416b707SBarry Smith PetscOptionItem amsopt; 101044ef3d73SBarry Smith PetscBool lset; 101153acd3b1SBarry Smith 101253acd3b1SBarry Smith PetscFunctionBegin; 1013e55864a3SBarry Smith if (!PetscOptionsObject->count) { 10149566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_REAL, &amsopt)); 10159566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscReal), &amsopt->data)); 1016a297a907SKarl Rupp 10170fdccdaeSBarry Smith *(PetscReal *)amsopt->data = currentvalue; 1018538aa990SBarry Smith } 10199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, &lset)); 102044ef3d73SBarry Smith if (set) *set = lset; 10211a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10229566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%g : %g>: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, lset && value ? (double)*value : (double)currentvalue, (double)currentvalue, text, ManSection(man))); 102353acd3b1SBarry Smith } 10243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102553acd3b1SBarry Smith } 102653acd3b1SBarry Smith 102788aa4217SBarry Smith /*MC 1028811af0c4SBarry Smith PetscOptionsScalar - Gets the `PetscScalar` value for a particular option in the database. 102953acd3b1SBarry Smith 1030811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 103153acd3b1SBarry Smith 103288aa4217SBarry Smith Synopsis: 103388aa4217SBarry Smith #include "petscsys.h" 10343a89f35bSSatish Balay PetscErrorCode PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 103588aa4217SBarry Smith 103653acd3b1SBarry Smith Input Parameters: 103753acd3b1SBarry Smith + opt - option name 103853acd3b1SBarry Smith . text - short string that describes the option 103953acd3b1SBarry Smith . man - manual page with additional information on option 10400fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 104188aa4217SBarry Smith $ PetscOptionsScalar(..., obj->value,&obj->value,...) or 10420fdccdaeSBarry Smith $ value = defaultvalue 10430fdccdaeSBarry Smith $ PetscOptionsScalar(..., value,&value,&flg); 10440fdccdaeSBarry Smith $ if (flg) { 10450fdccdaeSBarry Smith 1046d8d19677SJose E. Roman Output Parameters: 104753acd3b1SBarry Smith + value - the value to return 1048811af0c4SBarry Smith - flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 104953acd3b1SBarry Smith 105095452b02SPatrick Sanan Notes: 105195452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 10522efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 10532efd9cb1SBarry Smith 1054989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 1055989712b9SBarry Smith 1056811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 105753acd3b1SBarry Smith 10585a856986SBarry Smith Level: beginner 10595a856986SBarry Smith 1060db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1061db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1062db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1063db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1064c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1065db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1066db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 106788aa4217SBarry Smith M*/ 106888aa4217SBarry Smith 1069d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscScalar currentvalue, PetscScalar *value, PetscBool *set) 1070d71ae5a4SJacob Faibussowitsch { 107153acd3b1SBarry Smith PetscFunctionBegin; 107253acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 10739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal(opt, text, man, currentvalue, value, set)); 107453acd3b1SBarry Smith #else 10759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, set)); 107653acd3b1SBarry Smith #endif 10773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107853acd3b1SBarry Smith } 107953acd3b1SBarry Smith 108088aa4217SBarry Smith /*MC 108190d69ab7SBarry Smith PetscOptionsName - Determines if a particular option has been set in the database. This returns true whether the option is a number, string or boolean, even 108290d69ab7SBarry Smith its value is set to false. 108353acd3b1SBarry Smith 1084811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 108553acd3b1SBarry Smith 108688aa4217SBarry Smith Synopsis: 108788aa4217SBarry Smith #include "petscsys.h" 10883a89f35bSSatish Balay PetscErrorCode PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool *flg) 108988aa4217SBarry Smith 109053acd3b1SBarry Smith Input Parameters: 109153acd3b1SBarry Smith + opt - option name 109253acd3b1SBarry Smith . text - short string that describes the option 109353acd3b1SBarry Smith - man - manual page with additional information on option 109453acd3b1SBarry Smith 109553acd3b1SBarry Smith Output Parameter: 1096811af0c4SBarry Smith . flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 109753acd3b1SBarry Smith 109853acd3b1SBarry Smith Level: beginner 109953acd3b1SBarry Smith 1100811af0c4SBarry Smith Note: 1101811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 110253acd3b1SBarry Smith 1103db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1104db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1105db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1106db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1107c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1108db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1109db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 111088aa4217SBarry Smith M*/ 111188aa4217SBarry Smith 1112d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool *flg) 1113d71ae5a4SJacob Faibussowitsch { 11144416b707SBarry Smith PetscOptionItem amsopt; 111553acd3b1SBarry Smith 111653acd3b1SBarry Smith PetscFunctionBegin; 1117e55864a3SBarry Smith if (!PetscOptionsObject->count) { 11189566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL, &amsopt)); 11199566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool), &amsopt->data)); 1120a297a907SKarl Rupp 1121ace3abfcSBarry Smith *(PetscBool *)amsopt->data = PETSC_FALSE; 11221ae3d29cSBarry Smith } 11239566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, flg)); 1124e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 11259566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, text, ManSection(man))); 112653acd3b1SBarry Smith } 11273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112853acd3b1SBarry Smith } 112953acd3b1SBarry Smith 113088aa4217SBarry Smith /*MC 1131a264d7a6SBarry Smith PetscOptionsFList - Puts a list of option values that a single one may be selected from 113253acd3b1SBarry Smith 1133811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 113453acd3b1SBarry Smith 113588aa4217SBarry Smith Synopsis: 113688aa4217SBarry Smith #include "petscsys.h" 11373a89f35bSSatish Balay PetscErrorCode PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool *set) 113888aa4217SBarry Smith 113953acd3b1SBarry Smith Input Parameters: 114053acd3b1SBarry Smith + opt - option name 114153acd3b1SBarry Smith . text - short string that describes the option 114253acd3b1SBarry Smith . man - manual page with additional information on option 114353acd3b1SBarry Smith . list - the possible choices 11440fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 11450fdccdaeSBarry Smith $ PetscOptionsFlist(..., obj->value,value,len,&flg); 11460fdccdaeSBarry Smith $ if (flg) { 11473cc1e11dSBarry Smith - len - the length of the character array value 114853acd3b1SBarry Smith 1149d8d19677SJose E. Roman Output Parameters: 115053acd3b1SBarry Smith + value - the value to return 1151811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 115253acd3b1SBarry Smith 115353acd3b1SBarry Smith Level: intermediate 115453acd3b1SBarry Smith 115595452b02SPatrick Sanan Notes: 1156811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 115753acd3b1SBarry Smith 11582efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 11592efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 11602efd9cb1SBarry Smith 1161989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 1162989712b9SBarry Smith 1163811af0c4SBarry Smith See `PetscOptionsEList()` for when the choices are given in a string array 116453acd3b1SBarry Smith 116553acd3b1SBarry Smith To get a listing of all currently specified options, 1166811af0c4SBarry Smith see `PetscOptionsView()` or `PetscOptionsGetAll()` 116753acd3b1SBarry Smith 1168811af0c4SBarry Smith Developer Note: 1169811af0c4SBarry Smith This cannot check for invalid selection because of things like `MATAIJ` that are not included in the list 1170eabe10d7SBarry Smith 1171db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1172db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1173db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1174c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1175db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1176db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 117788aa4217SBarry Smith M*/ 117888aa4217SBarry Smith 1179d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsFList_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char ltext[], const char man[], PetscFunctionList list, const char currentvalue[], char value[], size_t len, PetscBool *set) 1180d71ae5a4SJacob Faibussowitsch { 11814416b707SBarry Smith PetscOptionItem amsopt; 118244ef3d73SBarry Smith PetscBool lset; 118353acd3b1SBarry Smith 118453acd3b1SBarry Smith PetscFunctionBegin; 11851a1499c8SBarry Smith if (!PetscOptionsObject->count) { 11869566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, ltext, man, OPTION_FLIST, &amsopt)); 118764facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 11889566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "", (char **)&amsopt->data)); 11893cc1e11dSBarry Smith amsopt->flist = list; 11903cc1e11dSBarry Smith } 11919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, len, &lset)); 119244ef3d73SBarry Smith if (set) *set = lset; 11931a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 11949566063dSJacob Faibussowitsch PetscCall(PetscFunctionListPrintTypes(PetscOptionsObject->comm, stdout, PetscOptionsObject->prefix, opt, ltext, man, list, currentvalue, lset && value ? value : currentvalue)); 119553acd3b1SBarry Smith } 11963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119753acd3b1SBarry Smith } 119853acd3b1SBarry Smith 119988aa4217SBarry Smith /*MC 120053acd3b1SBarry Smith PetscOptionsEList - Puts a list of option values that a single one may be selected from 120153acd3b1SBarry Smith 1202811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 120353acd3b1SBarry Smith 120488aa4217SBarry Smith Synopsis: 120588aa4217SBarry Smith #include "petscsys.h" 12063a89f35bSSatish Balay PetscErrorCode PetscOptionsEList(const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool *set) 120788aa4217SBarry Smith 120853acd3b1SBarry Smith Input Parameters: 120953acd3b1SBarry Smith + opt - option name 121053acd3b1SBarry Smith . ltext - short string that describes the option 121153acd3b1SBarry Smith . man - manual page with additional information on option 1212a264d7a6SBarry Smith . list - the possible choices (one of these must be selected, anything else is invalid) 121353acd3b1SBarry Smith . ntext - number of choices 12140fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 12150fdccdaeSBarry Smith $ PetscOptionsElist(..., obj->value,&value,&flg); 12160fdccdaeSBarry Smith $ if (flg) { 12170fdccdaeSBarry Smith 1218d8d19677SJose E. Roman Output Parameters: 121953acd3b1SBarry Smith + value - the index of the value to return 1220811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 122153acd3b1SBarry Smith 122253acd3b1SBarry Smith Level: intermediate 122353acd3b1SBarry Smith 122495452b02SPatrick Sanan Notes: 1225811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 122653acd3b1SBarry Smith 12272efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 12282efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 12292efd9cb1SBarry Smith 1230811af0c4SBarry Smith See `PetscOptionsFList()` for when the choices are given in a `PetscFunctionList()` 123153acd3b1SBarry Smith 1232db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1233db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1234db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1235c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1236db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1237db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEnum()` 123888aa4217SBarry Smith M*/ 123988aa4217SBarry Smith 1240d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsEList_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char ltext[], const char man[], const char *const *list, PetscInt ntext, const char currentvalue[], PetscInt *value, PetscBool *set) 1241d71ae5a4SJacob Faibussowitsch { 124253acd3b1SBarry Smith PetscInt i; 12434416b707SBarry Smith PetscOptionItem amsopt; 124444ef3d73SBarry Smith PetscBool lset; 124553acd3b1SBarry Smith 124653acd3b1SBarry Smith PetscFunctionBegin; 12471a1499c8SBarry Smith if (!PetscOptionsObject->count) { 12489566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, ltext, man, OPTION_ELIST, &amsopt)); 124964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 12509566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "", (char **)&amsopt->data)); 12519566063dSJacob Faibussowitsch PetscCall(PetscStrNArrayallocpy(ntext, list, (char ***)&amsopt->list)); 12521ae3d29cSBarry Smith amsopt->nlist = ntext; 12531ae3d29cSBarry Smith } 12549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEList(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, list, ntext, value, &lset)); 125544ef3d73SBarry Smith if (set) *set = lset; 12561a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 12579566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <now %s : formerly %s> %s (choose one of)", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, lset && value ? list[*value] : currentvalue, currentvalue, ltext)); 125848a46eb9SPierre Jolivet for (i = 0; i < ntext; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " %s", list[i])); 12599566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " (%s)\n", ManSection(man))); 126053acd3b1SBarry Smith } 12613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126253acd3b1SBarry Smith } 126353acd3b1SBarry Smith 126488aa4217SBarry Smith /*MC 1265acfcf0e5SJed Brown PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 1266d5649816SBarry Smith which at most a single value can be true. 126753acd3b1SBarry Smith 1268811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 126953acd3b1SBarry Smith 127088aa4217SBarry Smith Synopsis: 127188aa4217SBarry Smith #include "petscsys.h" 12723a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool *flg) 127388aa4217SBarry Smith 127453acd3b1SBarry Smith Input Parameters: 127553acd3b1SBarry Smith + opt - option name 127653acd3b1SBarry Smith . text - short string that describes the option 127753acd3b1SBarry Smith - man - manual page with additional information on option 127853acd3b1SBarry Smith 127953acd3b1SBarry Smith Output Parameter: 128053acd3b1SBarry Smith . flg - whether that option was set or not 128153acd3b1SBarry Smith 128253acd3b1SBarry Smith Level: intermediate 128353acd3b1SBarry Smith 128495452b02SPatrick Sanan Notes: 1285811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 128653acd3b1SBarry Smith 1287811af0c4SBarry Smith Must be followed by 0 or more P`etscOptionsBoolGroup()`s and `PetscOptionsBoolGroupEnd()` 128853acd3b1SBarry Smith 1289db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1290db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1291db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1292c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1293db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1294db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 129588aa4217SBarry Smith M*/ 129688aa4217SBarry Smith 1297d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool *flg) 1298d71ae5a4SJacob Faibussowitsch { 12994416b707SBarry Smith PetscOptionItem amsopt; 130053acd3b1SBarry Smith 130153acd3b1SBarry Smith PetscFunctionBegin; 1302e55864a3SBarry Smith if (!PetscOptionsObject->count) { 13039566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL, &amsopt)); 13049566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool), &amsopt->data)); 1305a297a907SKarl Rupp 1306ace3abfcSBarry Smith *(PetscBool *)amsopt->data = PETSC_FALSE; 13071ae3d29cSBarry Smith } 130868b16fdaSBarry Smith *flg = PETSC_FALSE; 13099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, flg, NULL)); 1310e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 13119566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " Pick at most one of -------------\n")); 13129566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, text, ManSection(man))); 131353acd3b1SBarry Smith } 13143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 131553acd3b1SBarry Smith } 131653acd3b1SBarry Smith 131788aa4217SBarry Smith /*MC 1318acfcf0e5SJed Brown PetscOptionsBoolGroup - One in a series of logical queries on the options database for 1319d5649816SBarry Smith which at most a single value can be true. 132053acd3b1SBarry Smith 1321811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 132253acd3b1SBarry Smith 132388aa4217SBarry Smith Synopsis: 132488aa4217SBarry Smith #include "petscsys.h" 13253a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool *flg) 132688aa4217SBarry Smith 132753acd3b1SBarry Smith Input Parameters: 132853acd3b1SBarry Smith + opt - option name 132953acd3b1SBarry Smith . text - short string that describes the option 133053acd3b1SBarry Smith - man - manual page with additional information on option 133153acd3b1SBarry Smith 133253acd3b1SBarry Smith Output Parameter: 1333811af0c4SBarry Smith . flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 133453acd3b1SBarry Smith 133553acd3b1SBarry Smith Level: intermediate 133653acd3b1SBarry Smith 133795452b02SPatrick Sanan Notes: 1338811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 133953acd3b1SBarry Smith 1340811af0c4SBarry Smith Must follow a `PetscOptionsBoolGroupBegin()` and preceded a `PetscOptionsBoolGroupEnd()` 134153acd3b1SBarry Smith 1342db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1343db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1344db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1345c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1346db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1347db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 134888aa4217SBarry Smith M*/ 134988aa4217SBarry Smith 1350d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool *flg) 1351d71ae5a4SJacob Faibussowitsch { 13524416b707SBarry Smith PetscOptionItem amsopt; 135353acd3b1SBarry Smith 135453acd3b1SBarry Smith PetscFunctionBegin; 1355e55864a3SBarry Smith if (!PetscOptionsObject->count) { 13569566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL, &amsopt)); 13579566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool), &amsopt->data)); 1358a297a907SKarl Rupp 1359ace3abfcSBarry Smith *(PetscBool *)amsopt->data = PETSC_FALSE; 13601ae3d29cSBarry Smith } 136117326d04SJed Brown *flg = PETSC_FALSE; 13629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, flg, NULL)); 1363e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 13649566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, text, ManSection(man))); 136553acd3b1SBarry Smith } 13663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136753acd3b1SBarry Smith } 136853acd3b1SBarry Smith 136988aa4217SBarry Smith /*MC 1370acfcf0e5SJed Brown PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 1371d5649816SBarry Smith which at most a single value can be true. 137253acd3b1SBarry Smith 1373811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 137453acd3b1SBarry Smith 137588aa4217SBarry Smith Synopsis: 137688aa4217SBarry Smith #include "petscsys.h" 13773a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool *flg) 137888aa4217SBarry Smith 137953acd3b1SBarry Smith Input Parameters: 138053acd3b1SBarry Smith + opt - option name 138153acd3b1SBarry Smith . text - short string that describes the option 138253acd3b1SBarry Smith - man - manual page with additional information on option 138353acd3b1SBarry Smith 138453acd3b1SBarry Smith Output Parameter: 1385811af0c4SBarry Smith . flg - `PETSC_TRUE` if found, else `PETSC_FALSE` 138653acd3b1SBarry Smith 138753acd3b1SBarry Smith Level: intermediate 138853acd3b1SBarry Smith 138995452b02SPatrick Sanan Notes: 1390811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 139153acd3b1SBarry Smith 1392811af0c4SBarry Smith Must follow a `PetscOptionsBoolGroupBegin()` 139353acd3b1SBarry Smith 1394db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1395db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1396db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1397c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1398db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1399db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 140088aa4217SBarry Smith M*/ 140188aa4217SBarry Smith 1402d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool *flg) 1403d71ae5a4SJacob Faibussowitsch { 14044416b707SBarry Smith PetscOptionItem amsopt; 140553acd3b1SBarry Smith 140653acd3b1SBarry Smith PetscFunctionBegin; 1407e55864a3SBarry Smith if (!PetscOptionsObject->count) { 14089566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL, &amsopt)); 14099566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool), &amsopt->data)); 1410a297a907SKarl Rupp 1411ace3abfcSBarry Smith *(PetscBool *)amsopt->data = PETSC_FALSE; 14121ae3d29cSBarry Smith } 141317326d04SJed Brown *flg = PETSC_FALSE; 14149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, flg, NULL)); 1415e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 14169566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, text, ManSection(man))); 141753acd3b1SBarry Smith } 14183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141953acd3b1SBarry Smith } 142053acd3b1SBarry Smith 142188aa4217SBarry Smith /*MC 1422acfcf0e5SJed Brown PetscOptionsBool - Determines if a particular option is in the database with a true or false 142353acd3b1SBarry Smith 1424811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 142553acd3b1SBarry Smith 142688aa4217SBarry Smith Synopsis: 142788aa4217SBarry Smith #include "petscsys.h" 14283a89f35bSSatish Balay PetscErrorCode PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 142988aa4217SBarry Smith 143053acd3b1SBarry Smith Input Parameters: 143153acd3b1SBarry Smith + opt - option name 143253acd3b1SBarry Smith . text - short string that describes the option 1433868c398cSBarry Smith . man - manual page with additional information on option 143494ae4db5SBarry Smith - currentvalue - the current value 143553acd3b1SBarry Smith 1436d8d19677SJose E. Roman Output Parameters: 1437811af0c4SBarry Smith + flg -` PETSC_TRUE` or `PETSC_FALSE` 1438811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 143953acd3b1SBarry Smith 14402efd9cb1SBarry Smith Notes: 1441811af0c4SBarry Smith TRUE, true, YES, yes, nostring, and 1 all translate to `PETSC_TRUE` 1442811af0c4SBarry Smith FALSE, false, NO, no, and 0 all translate to `PETSC_FALSE` 14432efd9cb1SBarry Smith 1444811af0c4SBarry Smith If the option is given, but no value is provided, then flg and set are both given the value `PETSC_TRUE`. That is -requested_bool 14452efd9cb1SBarry Smith is equivalent to -requested_bool true 14462efd9cb1SBarry Smith 14472efd9cb1SBarry Smith If the user does not supply the option at all flg is NOT changed. Thus 14482efd9cb1SBarry Smith you should ALWAYS initialize the flg if you access it without first checking if the set flag is true. 14492efd9cb1SBarry Smith 1450811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 145153acd3b1SBarry Smith 14525a856986SBarry Smith Level: beginner 14535a856986SBarry Smith 1454db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1455db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1456db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 1457db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1458c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1459db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1460db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 146188aa4217SBarry Smith M*/ 146288aa4217SBarry Smith 1463d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool currentvalue, PetscBool *flg, PetscBool *set) 1464d71ae5a4SJacob Faibussowitsch { 1465ace3abfcSBarry Smith PetscBool iset; 14664416b707SBarry Smith PetscOptionItem amsopt; 146753acd3b1SBarry Smith 146853acd3b1SBarry Smith PetscFunctionBegin; 1469e55864a3SBarry Smith if (!PetscOptionsObject->count) { 14709566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL, &amsopt)); 14719566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool), &amsopt->data)); 1472a297a907SKarl Rupp 147394ae4db5SBarry Smith *(PetscBool *)amsopt->data = currentvalue; 1474af6d86caSBarry Smith } 14759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, flg, &iset)); 147653acd3b1SBarry Smith if (set) *set = iset; 14771a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 147844ef3d73SBarry Smith const char *v = PetscBools[currentvalue], *vn = PetscBools[iset && flg ? *flg : currentvalue]; 14799566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s: <%s : %s> %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, v, vn, text, ManSection(man))); 148053acd3b1SBarry Smith } 14813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148253acd3b1SBarry Smith } 148353acd3b1SBarry Smith 148488aa4217SBarry Smith /*MC 148553acd3b1SBarry Smith PetscOptionsRealArray - Gets an array of double values for a particular 148653acd3b1SBarry Smith option in the database. The values must be separated with commas with 148753acd3b1SBarry Smith no intervening spaces. 148853acd3b1SBarry Smith 1489811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 149053acd3b1SBarry Smith 149188aa4217SBarry Smith Synopsis: 149288aa4217SBarry Smith #include "petscsys.h" 14933a89f35bSSatish Balay PetscErrorCode PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 149488aa4217SBarry Smith 149553acd3b1SBarry Smith Input Parameters: 149653acd3b1SBarry Smith + opt - the option one is seeking 149753acd3b1SBarry Smith . text - short string describing option 149853acd3b1SBarry Smith . man - manual page for option 1499c89788bdSBarry Smith - n - maximum number of values that value has room for 150053acd3b1SBarry Smith 1501d8d19677SJose E. Roman Output Parameters: 150253acd3b1SBarry Smith + value - location to copy values 1503c89788bdSBarry Smith . n - actual number of values found 1504811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 150553acd3b1SBarry Smith 150653acd3b1SBarry Smith Level: beginner 150753acd3b1SBarry Smith 1508811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 150953acd3b1SBarry Smith 1510db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1511db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1512db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1513c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1514db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1515db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 151688aa4217SBarry Smith M*/ 151788aa4217SBarry Smith 1518d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscReal value[], PetscInt *n, PetscBool *set) 1519d71ae5a4SJacob Faibussowitsch { 152053acd3b1SBarry Smith PetscInt i; 15214416b707SBarry Smith PetscOptionItem amsopt; 152253acd3b1SBarry Smith 152353acd3b1SBarry Smith PetscFunctionBegin; 1524e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1525e26ddf31SBarry Smith PetscReal *vals; 1526e26ddf31SBarry Smith 15279566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_REAL_ARRAY, &amsopt)); 15289566063dSJacob Faibussowitsch PetscCall(PetscMalloc((*n) * sizeof(PetscReal), &amsopt->data)); 1529e26ddf31SBarry Smith vals = (PetscReal *)amsopt->data; 1530e26ddf31SBarry Smith for (i = 0; i < *n; i++) vals[i] = value[i]; 1531e26ddf31SBarry Smith amsopt->arraylength = *n; 1532e26ddf31SBarry Smith } 15339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, n, set)); 1534e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 15359566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%g", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, (double)value[0])); 153648a46eb9SPierre Jolivet for (i = 1; i < *n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ",%g", (double)value[i])); 15379566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ">: %s (%s)\n", text, ManSection(man))); 153853acd3b1SBarry Smith } 15393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 154053acd3b1SBarry Smith } 154153acd3b1SBarry Smith 154288aa4217SBarry Smith /*MC 1543811af0c4SBarry Smith PetscOptionsScalarArray - Gets an array of `PetscScalar` values for a particular 1544050cccc3SHong Zhang option in the database. The values must be separated with commas with 1545050cccc3SHong Zhang no intervening spaces. 1546050cccc3SHong Zhang 1547811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 1548050cccc3SHong Zhang 154988aa4217SBarry Smith Synopsis: 155088aa4217SBarry Smith #include "petscsys.h" 15513a89f35bSSatish Balay PetscErrorCode PetscOptionsScalarArray(const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool *set) 155288aa4217SBarry Smith 1553050cccc3SHong Zhang Input Parameters: 1554050cccc3SHong Zhang + opt - the option one is seeking 1555050cccc3SHong Zhang . text - short string describing option 1556050cccc3SHong Zhang . man - manual page for option 1557c89788bdSBarry Smith - n - maximum number of values allowed in the value array 1558050cccc3SHong Zhang 1559d8d19677SJose E. Roman Output Parameters: 1560050cccc3SHong Zhang + value - location to copy values 1561c89788bdSBarry Smith . n - actual number of values found 1562811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 1563050cccc3SHong Zhang 1564050cccc3SHong Zhang Level: beginner 1565050cccc3SHong Zhang 1566811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 1567050cccc3SHong Zhang 1568db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1569db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1570db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1571c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1572db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1573db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 157488aa4217SBarry Smith M*/ 157588aa4217SBarry Smith 1576d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscScalar value[], PetscInt *n, PetscBool *set) 1577d71ae5a4SJacob Faibussowitsch { 1578050cccc3SHong Zhang PetscInt i; 15794416b707SBarry Smith PetscOptionItem amsopt; 1580050cccc3SHong Zhang 1581050cccc3SHong Zhang PetscFunctionBegin; 1582050cccc3SHong Zhang if (!PetscOptionsObject->count) { 1583050cccc3SHong Zhang PetscScalar *vals; 1584050cccc3SHong Zhang 15859566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_SCALAR_ARRAY, &amsopt)); 15869566063dSJacob Faibussowitsch PetscCall(PetscMalloc((*n) * sizeof(PetscScalar), &amsopt->data)); 1587050cccc3SHong Zhang vals = (PetscScalar *)amsopt->data; 1588050cccc3SHong Zhang for (i = 0; i < *n; i++) vals[i] = value[i]; 1589050cccc3SHong Zhang amsopt->arraylength = *n; 1590050cccc3SHong Zhang } 15919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalarArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, n, set)); 1592050cccc3SHong Zhang if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 15939566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%g+%gi", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, (double)PetscRealPart(value[0]), (double)PetscImaginaryPart(value[0]))); 159448a46eb9SPierre Jolivet for (i = 1; i < *n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ",%g+%gi", (double)PetscRealPart(value[i]), (double)PetscImaginaryPart(value[i]))); 15959566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ">: %s (%s)\n", text, ManSection(man))); 1596050cccc3SHong Zhang } 15973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1598050cccc3SHong Zhang } 159953acd3b1SBarry Smith 160088aa4217SBarry Smith /*MC 160153acd3b1SBarry Smith PetscOptionsIntArray - Gets an array of integers for a particular 1602b32a342fSShri Abhyankar option in the database. 160353acd3b1SBarry Smith 1604811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 160553acd3b1SBarry Smith 160688aa4217SBarry Smith Synopsis: 160788aa4217SBarry Smith #include "petscsys.h" 16083a89f35bSSatish Balay PetscErrorCode PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 160988aa4217SBarry Smith 161053acd3b1SBarry Smith Input Parameters: 161153acd3b1SBarry Smith + opt - the option one is seeking 161253acd3b1SBarry Smith . text - short string describing option 161353acd3b1SBarry Smith . man - manual page for option 1614f8a50e2bSBarry Smith - n - maximum number of values 161553acd3b1SBarry Smith 1616d8d19677SJose E. Roman Output Parameters: 161753acd3b1SBarry Smith + value - location to copy values 1618f8a50e2bSBarry Smith . n - actual number of values found 1619811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 162053acd3b1SBarry Smith 162153acd3b1SBarry Smith Level: beginner 162253acd3b1SBarry Smith 162353acd3b1SBarry Smith Notes: 1624b32a342fSShri Abhyankar The array can be passed as 1625811af0c4SBarry Smith + a comma separated list - 0,1,2,3,4,5,6,7 1626811af0c4SBarry Smith . a range (start\-end+1) - 0-8 1627811af0c4SBarry Smith . a range with given increment (start\-end+1:inc) - 0-7:2 1628811af0c4SBarry Smith - a combination of values and ranges separated by commas - 0,1-8,8-15:2 1629b32a342fSShri Abhyankar 1630b32a342fSShri Abhyankar There must be no intervening spaces between the values. 163153acd3b1SBarry Smith 1632811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 163353acd3b1SBarry Smith 1634db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1635db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1636db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1637c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1638db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1639db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRealArray()` 164088aa4217SBarry Smith M*/ 164188aa4217SBarry Smith 1642d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscInt value[], PetscInt *n, PetscBool *set) 1643d71ae5a4SJacob Faibussowitsch { 164453acd3b1SBarry Smith PetscInt i; 16454416b707SBarry Smith PetscOptionItem amsopt; 164653acd3b1SBarry Smith 164753acd3b1SBarry Smith PetscFunctionBegin; 1648e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1649e26ddf31SBarry Smith PetscInt *vals; 1650e26ddf31SBarry Smith 16519566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_INT_ARRAY, &amsopt)); 16529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n, (PetscInt **)&amsopt->data)); 1653e26ddf31SBarry Smith vals = (PetscInt *)amsopt->data; 1654e26ddf31SBarry Smith for (i = 0; i < *n; i++) vals[i] = value[i]; 1655e26ddf31SBarry Smith amsopt->arraylength = *n; 1656e26ddf31SBarry Smith } 16579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, n, set)); 1658e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 16599566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%" PetscInt_FMT, PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, value[0])); 166048a46eb9SPierre Jolivet for (i = 1; i < *n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ",%" PetscInt_FMT, value[i])); 16619566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ">: %s (%s)\n", text, ManSection(man))); 166253acd3b1SBarry Smith } 16633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166453acd3b1SBarry Smith } 166553acd3b1SBarry Smith 166688aa4217SBarry Smith /*MC 166753acd3b1SBarry Smith PetscOptionsStringArray - Gets an array of string values for a particular 166853acd3b1SBarry Smith option in the database. The values must be separated with commas with 166953acd3b1SBarry Smith no intervening spaces. 167053acd3b1SBarry Smith 1671cf53795eSBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()`; No Fortran Support 167253acd3b1SBarry Smith 167388aa4217SBarry Smith Synopsis: 167488aa4217SBarry Smith #include "petscsys.h" 16753a89f35bSSatish Balay PetscErrorCode PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 167688aa4217SBarry Smith 167753acd3b1SBarry Smith Input Parameters: 167853acd3b1SBarry Smith + opt - the option one is seeking 167953acd3b1SBarry Smith . text - short string describing option 168053acd3b1SBarry Smith . man - manual page for option 168153acd3b1SBarry Smith - nmax - maximum number of strings 168253acd3b1SBarry Smith 1683d8d19677SJose E. Roman Output Parameters: 168453acd3b1SBarry Smith + value - location to copy strings 168553acd3b1SBarry Smith . nmax - actual number of strings found 1686811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 168753acd3b1SBarry Smith 168853acd3b1SBarry Smith Level: beginner 168953acd3b1SBarry Smith 169053acd3b1SBarry Smith Notes: 169153acd3b1SBarry Smith The user should pass in an array of pointers to char, to hold all the 169253acd3b1SBarry Smith strings returned by this function. 169353acd3b1SBarry Smith 169453acd3b1SBarry Smith The user is responsible for deallocating the strings that are 1695cf53795eSBarry Smith returned. 169653acd3b1SBarry Smith 1697811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 169853acd3b1SBarry Smith 1699db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1700db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1701db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1702c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1703db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1704db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 170588aa4217SBarry Smith M*/ 170688aa4217SBarry Smith 1707d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], char *value[], PetscInt *nmax, PetscBool *set) 1708d71ae5a4SJacob Faibussowitsch { 17094416b707SBarry Smith PetscOptionItem amsopt; 171053acd3b1SBarry Smith 171153acd3b1SBarry Smith PetscFunctionBegin; 1712e55864a3SBarry Smith if (!PetscOptionsObject->count) { 17139566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_STRING_ARRAY, &amsopt)); 17149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nmax, (char **)&amsopt->data)); 1715a297a907SKarl Rupp 17161ae3d29cSBarry Smith amsopt->arraylength = *nmax; 17171ae3d29cSBarry Smith } 17189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, nmax, set)); 1719e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 17209566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <string1,string2,...>: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, text, ManSection(man))); 172153acd3b1SBarry Smith } 17223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172353acd3b1SBarry Smith } 172453acd3b1SBarry Smith 172588aa4217SBarry Smith /*MC 1726acfcf0e5SJed Brown PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 1727e2446a98SMatthew Knepley option in the database. The values must be separated with commas with 1728e2446a98SMatthew Knepley no intervening spaces. 1729e2446a98SMatthew Knepley 1730811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 1731e2446a98SMatthew Knepley 173288aa4217SBarry Smith Synopsis: 173388aa4217SBarry Smith #include "petscsys.h" 17343a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 173588aa4217SBarry Smith 1736e2446a98SMatthew Knepley Input Parameters: 1737e2446a98SMatthew Knepley + opt - the option one is seeking 1738e2446a98SMatthew Knepley . text - short string describing option 1739e2446a98SMatthew Knepley . man - manual page for option 1740c89788bdSBarry Smith - n - maximum number of values allowed in the value array 1741e2446a98SMatthew Knepley 1742d8d19677SJose E. Roman Output Parameters: 1743e2446a98SMatthew Knepley + value - location to copy values 1744c89788bdSBarry Smith . n - actual number of values found 1745811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 1746e2446a98SMatthew Knepley 1747e2446a98SMatthew Knepley Level: beginner 1748e2446a98SMatthew Knepley 1749e2446a98SMatthew Knepley Notes: 1750e2446a98SMatthew Knepley The user should pass in an array of doubles 1751e2446a98SMatthew Knepley 1752811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 1753e2446a98SMatthew Knepley 1754db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1755db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1756db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1757c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1758db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1759db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 176088aa4217SBarry Smith M*/ 176188aa4217SBarry Smith 1762d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscBool value[], PetscInt *n, PetscBool *set) 1763d71ae5a4SJacob Faibussowitsch { 1764e2446a98SMatthew Knepley PetscInt i; 17654416b707SBarry Smith PetscOptionItem amsopt; 1766e2446a98SMatthew Knepley 1767e2446a98SMatthew Knepley PetscFunctionBegin; 1768e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1769ace3abfcSBarry Smith PetscBool *vals; 17701ae3d29cSBarry Smith 17719566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_BOOL_ARRAY, &amsopt)); 17729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n, (PetscBool **)&amsopt->data)); 1773ace3abfcSBarry Smith vals = (PetscBool *)amsopt->data; 17741ae3d29cSBarry Smith for (i = 0; i < *n; i++) vals[i] = value[i]; 17751ae3d29cSBarry Smith amsopt->arraylength = *n; 17761ae3d29cSBarry Smith } 17779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBoolArray(PetscOptionsObject->options, PetscOptionsObject->prefix, opt, value, n, set)); 1778e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 17799566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%d", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, value[0])); 178048a46eb9SPierre Jolivet for (i = 1; i < *n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ",%d", value[i])); 17819566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, ">: %s (%s)\n", text, ManSection(man))); 1782e2446a98SMatthew Knepley } 17833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1784e2446a98SMatthew Knepley } 1785e2446a98SMatthew Knepley 178688aa4217SBarry Smith /*MC 1787d1da0b69SBarry Smith PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user 17888cc676e6SMatthew G Knepley 1789811af0c4SBarry Smith Logically Collective on the communicator passed in `PetscOptionsBegin()` 17908cc676e6SMatthew G Knepley 179188aa4217SBarry Smith Synopsis: 179288aa4217SBarry Smith #include "petscsys.h" 17933a89f35bSSatish Balay PetscErrorCode PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 179488aa4217SBarry Smith 17958cc676e6SMatthew G Knepley Input Parameters: 17968cc676e6SMatthew G Knepley + opt - option name 17978cc676e6SMatthew G Knepley . text - short string that describes the option 17988cc676e6SMatthew G Knepley - man - manual page with additional information on option 17998cc676e6SMatthew G Knepley 1800d8d19677SJose E. Roman Output Parameters: 18018cc676e6SMatthew G Knepley + viewer - the viewer 18027962402dSFande Kong . format - the PetscViewerFormat requested by the user, pass NULL if not needed 1803811af0c4SBarry Smith - set - `PETSC_TRUE` if found, else `PETSC_FALSE` 18048cc676e6SMatthew G Knepley 18058cc676e6SMatthew G Knepley Level: beginner 18068cc676e6SMatthew G Knepley 180795452b02SPatrick Sanan Notes: 1808811af0c4SBarry Smith Must be between a `PetscOptionsBegin()` and a `PetscOptionsEnd()` 18098cc676e6SMatthew G Knepley 1810811af0c4SBarry Smith See `PetscOptionsGetViewer()` for the format of the supplied viewer and its options 18118cc676e6SMatthew G Knepley 1812db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1813db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1814db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1815db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1816c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1817db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1818db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 181988aa4217SBarry Smith M*/ 182088aa4217SBarry Smith 1821d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject, const char opt[], const char text[], const char man[], PetscViewer *viewer, PetscViewerFormat *format, PetscBool *set) 1822d71ae5a4SJacob Faibussowitsch { 18234416b707SBarry Smith PetscOptionItem amsopt; 18248cc676e6SMatthew G Knepley 18258cc676e6SMatthew G Knepley PetscFunctionBegin; 18261a1499c8SBarry Smith if (!PetscOptionsObject->count) { 18279566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject, opt, text, man, OPTION_STRING, &amsopt)); 182864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 18299566063dSJacob Faibussowitsch PetscCall(PetscStrdup("", (char **)&amsopt->data)); 18308cc676e6SMatthew G Knepley } 18319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscOptionsObject->comm, PetscOptionsObject->options, PetscOptionsObject->prefix, opt, viewer, format, set)); 1832e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 18339566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm, " -%s%s <%s>: %s (%s)\n", PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "", opt + 1, "", text, ManSection(man))); 18348cc676e6SMatthew G Knepley } 18353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18368cc676e6SMatthew G Knepley } 1837