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 */ 224416b707SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[]) 2353acd3b1SBarry Smith { 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) { 41e55864a3SBarry Smith if (!PetscOptionsObject->alreadyprinted) { 429566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(comm,"----------------------------------------\n%s:\n",title)); 4353acd3b1SBarry Smith } 4461b37b28SSatish Balay } 4553acd3b1SBarry Smith PetscFunctionReturn(0); 4653acd3b1SBarry Smith } 4753acd3b1SBarry Smith 483194b578SJed Brown /* 493194b578SJed Brown Handles setting up the data structure in a call to PetscObjectOptionsBegin() 503194b578SJed Brown */ 514416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,PetscObject obj) 523194b578SJed Brown { 533194b578SJed Brown char title[256]; 543194b578SJed Brown PetscBool flg; 553194b578SJed Brown 563194b578SJed Brown PetscFunctionBegin; 575f80ce2aSJacob Faibussowitsch PetscValidPointer(PetscOptionsObject,1); 58064a246eSJacob Faibussowitsch PetscValidHeader(obj,2); 59e55864a3SBarry Smith PetscOptionsObject->object = obj; 60e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = obj->optionsprinted; 61a297a907SKarl Rupp 629566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(obj->description,obj->class_name,&flg)); 639566063dSJacob Faibussowitsch if (flg) PetscCall(PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name)); 649566063dSJacob Faibussowitsch else PetscCall(PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec)); 663194b578SJed Brown PetscFunctionReturn(0); 673194b578SJed Brown } 683194b578SJed Brown 6953acd3b1SBarry Smith /* 7053acd3b1SBarry Smith Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd() 7153acd3b1SBarry Smith */ 724416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt) 7353acd3b1SBarry Smith { 744416b707SBarry Smith PetscOptionItem next; 753be6e4c3SJed Brown PetscBool valid; 7653acd3b1SBarry Smith 7753acd3b1SBarry Smith PetscFunctionBegin; 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsValidKey(opt,&valid)); 795f80ce2aSJacob Faibussowitsch PetscCheck(valid,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt); 803be6e4c3SJed Brown 819566063dSJacob Faibussowitsch PetscCall(PetscNew(amsopt)); 8202c9f0b5SLisandro Dalcin (*amsopt)->next = NULL; 8353acd3b1SBarry Smith (*amsopt)->set = PETSC_FALSE; 846356e834SBarry Smith (*amsopt)->type = t; 8502c9f0b5SLisandro Dalcin (*amsopt)->data = NULL; 8661b37b28SSatish Balay 879566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(text,&(*amsopt)->text)); 889566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(opt,&(*amsopt)->option)); 899566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(man,&(*amsopt)->man)); 9053acd3b1SBarry Smith 91e55864a3SBarry Smith if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt; 92a297a907SKarl Rupp else { 93e55864a3SBarry Smith next = PetscOptionsObject->next; 9453acd3b1SBarry Smith while (next->next) next = next->next; 9553acd3b1SBarry Smith next->next = *amsopt; 9653acd3b1SBarry Smith } 9753acd3b1SBarry Smith PetscFunctionReturn(0); 9853acd3b1SBarry Smith } 9953acd3b1SBarry Smith 100aee2cecaSBarry Smith /* 1013fc1eb6aSBarry Smith PetscScanString - Gets user input via stdin from process and broadcasts to all processes 1023fc1eb6aSBarry Smith 103d083f849SBarry Smith Collective 1043fc1eb6aSBarry Smith 1053fc1eb6aSBarry Smith Input Parameters: 1063fc1eb6aSBarry Smith + commm - communicator for the broadcast, must be PETSC_COMM_WORLD 1073fc1eb6aSBarry Smith . n - length of the string, must be the same on all processes 1083fc1eb6aSBarry Smith - str - location to store input 109aee2cecaSBarry Smith 110aee2cecaSBarry Smith Bugs: 111aee2cecaSBarry Smith . Assumes process 0 of the given communicator has access to stdin 112aee2cecaSBarry Smith 113aee2cecaSBarry Smith */ 1143fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm, size_t n, char str[]) 115aee2cecaSBarry Smith { 1163fc1eb6aSBarry Smith PetscMPIInt rank,nm; 117aee2cecaSBarry Smith 118aee2cecaSBarry Smith PetscFunctionBegin; 1199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 120dd400576SPatrick Sanan if (rank == 0) { 1215f80ce2aSJacob Faibussowitsch char c = (char)getchar(); 1225f80ce2aSJacob Faibussowitsch size_t i = 0; 1235f80ce2aSJacob Faibussowitsch 124aee2cecaSBarry Smith while (c != '\n' && i < n-1) { 125aee2cecaSBarry Smith str[i++] = c; 126aee2cecaSBarry Smith c = (char)getchar(); 127aee2cecaSBarry Smith } 128aee2cecaSBarry Smith str[i] = 0; 129aee2cecaSBarry Smith } 1309566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(n,&nm)); 1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(str,nm,MPI_CHAR,0,comm)); 132aee2cecaSBarry Smith PetscFunctionReturn(0); 133aee2cecaSBarry Smith } 134aee2cecaSBarry Smith 1355b02f95dSBarry Smith /* 1365b02f95dSBarry Smith This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy() 1375b02f95dSBarry Smith */ 1385b02f95dSBarry Smith static PetscErrorCode PetscStrdup(const char s[], char *t[]) 1395b02f95dSBarry Smith { 14002c9f0b5SLisandro Dalcin char *tmp = NULL; 1415b02f95dSBarry Smith 1425b02f95dSBarry Smith PetscFunctionBegin; 1435b02f95dSBarry Smith if (s) { 1445f80ce2aSJacob Faibussowitsch size_t len; 1455f80ce2aSJacob Faibussowitsch 1469566063dSJacob Faibussowitsch PetscCall(PetscStrlen(s,&len)); 1475f80ce2aSJacob Faibussowitsch tmp = (char*) malloc((len+1)*sizeof(*tmp)); 1485f80ce2aSJacob Faibussowitsch PetscCheck(tmp,PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string"); 1499566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(tmp,s)); 1505b02f95dSBarry Smith } 1515b02f95dSBarry Smith *t = tmp; 1525b02f95dSBarry Smith PetscFunctionReturn(0); 1535b02f95dSBarry Smith } 1545b02f95dSBarry Smith 155aee2cecaSBarry Smith /* 1563cc1e11dSBarry Smith PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime 157aee2cecaSBarry Smith 15895452b02SPatrick Sanan Notes: 15995452b02SPatrick Sanan this isn't really practical, it is just to demonstrate the principle 160aee2cecaSBarry Smith 1617781c08eSBarry Smith A carriage return indicates no change from the default; but this like -ksp_monitor <stdout> the default is actually not stdout the default 1627781c08eSBarry Smith is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug? 1637781c08eSBarry Smith 164aee2cecaSBarry Smith Bugs: 1657781c08eSBarry Smith + All processes must traverse through the exact same set of option queries due to the call to PetscScanString() 1663cc1e11dSBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 167aee2cecaSBarry Smith - Only works for PetscInt == int, PetscReal == double etc 168aee2cecaSBarry Smith 16995452b02SPatrick Sanan Developer Notes: 17095452b02SPatrick Sanan Normally the GUI that presents the options the user and retrieves the values would be running in a different 1713cc1e11dSBarry Smith address space and communicating with the PETSc program 1723cc1e11dSBarry Smith 173aee2cecaSBarry Smith */ 1744416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject) 1756356e834SBarry Smith { 1764416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 1776356e834SBarry Smith char str[512]; 1787781c08eSBarry Smith PetscBool bid; 179a4404d99SBarry Smith PetscReal ir,*valr; 180330cf3c9SBarry Smith PetscInt *vald; 181330cf3c9SBarry Smith size_t i; 1826356e834SBarry Smith 1832a409bb0SBarry Smith PetscFunctionBegin; 1849566063dSJacob Faibussowitsch PetscCall((*PetscPrintf)(PETSC_COMM_WORLD,"%s --------------------\n",PetscOptionsObject->title)); 1856356e834SBarry Smith while (next) { 1866356e834SBarry Smith switch (next->type) { 1876356e834SBarry Smith case OPTION_HEAD: 1886356e834SBarry Smith break; 189e26ddf31SBarry Smith case OPTION_INT_ARRAY: 1909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1)); 191e26ddf31SBarry Smith vald = (PetscInt*) next->data; 192e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 1939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT,vald[i])); 194e26ddf31SBarry Smith if (i < next->arraylength-1) { 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,",")); 196e26ddf31SBarry Smith } 197e26ddf31SBarry Smith } 1989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man)); 1999566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 200e26ddf31SBarry Smith if (str[0]) { 201e26ddf31SBarry Smith PetscToken token; 202e26ddf31SBarry Smith PetscInt n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end; 203e26ddf31SBarry Smith size_t len; 204e26ddf31SBarry Smith char *value; 205ace3abfcSBarry Smith PetscBool foundrange; 206e26ddf31SBarry Smith 207e26ddf31SBarry Smith next->set = PETSC_TRUE; 208e26ddf31SBarry Smith value = str; 2099566063dSJacob Faibussowitsch PetscCall(PetscTokenCreate(value,',',&token)); 2109566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token,&value)); 211e26ddf31SBarry Smith while (n < nmax) { 212e26ddf31SBarry Smith if (!value) break; 213e26ddf31SBarry Smith 214e26ddf31SBarry Smith /* look for form d-D where d and D are integers */ 215e26ddf31SBarry Smith foundrange = PETSC_FALSE; 2169566063dSJacob Faibussowitsch PetscCall(PetscStrlen(value,&len)); 217e26ddf31SBarry Smith if (value[0] == '-') i=2; 218e26ddf31SBarry Smith else i=1; 219330cf3c9SBarry Smith for (;i<len; i++) { 220e26ddf31SBarry Smith if (value[i] == '-') { 22108401ef6SPierre Jolivet PetscCheck(i != len-1,PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %" PetscInt_FMT "-th array entry %s",n,value); 222e26ddf31SBarry Smith value[i] = 0; 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value,&start)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value+i+1,&end)); 22508401ef6SPierre 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); 226cc73adaaSBarry 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); 227e26ddf31SBarry Smith for (; start<end; start++) { 228e26ddf31SBarry Smith *dvalue = start; dvalue++;n++; 229e26ddf31SBarry Smith } 230e26ddf31SBarry Smith foundrange = PETSC_TRUE; 231e26ddf31SBarry Smith break; 232e26ddf31SBarry Smith } 233e26ddf31SBarry Smith } 234e26ddf31SBarry Smith if (!foundrange) { 2359566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToInt(value,dvalue)); 236e26ddf31SBarry Smith dvalue++; 237e26ddf31SBarry Smith n++; 238e26ddf31SBarry Smith } 2399566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token,&value)); 240e26ddf31SBarry Smith } 2419566063dSJacob Faibussowitsch PetscCall(PetscTokenDestroy(&token)); 242e26ddf31SBarry Smith } 243e26ddf31SBarry Smith break; 244e26ddf31SBarry Smith case OPTION_REAL_ARRAY: 2459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1)); 246e26ddf31SBarry Smith valr = (PetscReal*) next->data; 247e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 2489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%g",(double)valr[i])); 249e26ddf31SBarry Smith if (i < next->arraylength-1) { 2509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,",")); 251e26ddf31SBarry Smith } 252e26ddf31SBarry Smith } 2539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man)); 2549566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 255e26ddf31SBarry Smith if (str[0]) { 256e26ddf31SBarry Smith PetscToken token; 257e26ddf31SBarry Smith PetscInt n = 0,nmax = next->arraylength; 258e26ddf31SBarry Smith PetscReal *dvalue = (PetscReal*)next->data; 259e26ddf31SBarry Smith char *value; 260e26ddf31SBarry Smith 261e26ddf31SBarry Smith next->set = PETSC_TRUE; 262e26ddf31SBarry Smith value = str; 2639566063dSJacob Faibussowitsch PetscCall(PetscTokenCreate(value,',',&token)); 2649566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token,&value)); 265e26ddf31SBarry Smith while (n < nmax) { 266e26ddf31SBarry Smith if (!value) break; 2679566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToReal(value,dvalue)); 268e26ddf31SBarry Smith dvalue++; 269e26ddf31SBarry Smith n++; 2709566063dSJacob Faibussowitsch PetscCall(PetscTokenFind(token,&value)); 271e26ddf31SBarry Smith } 2729566063dSJacob Faibussowitsch PetscCall(PetscTokenDestroy(&token)); 273e26ddf31SBarry Smith } 274e26ddf31SBarry Smith break; 2756356e834SBarry Smith case OPTION_INT: 2769566063dSJacob 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)); 2779566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 2783fc1eb6aSBarry Smith if (str[0]) { 279d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 280d25d7f95SJed Brown long long lid; 281d25d7f95SJed Brown sscanf(str,"%lld",&lid); 282cc73adaaSBarry 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); 283c272547aSJed Brown #else 284d25d7f95SJed Brown long lid; 285d25d7f95SJed Brown sscanf(str,"%ld",&lid); 286cc73adaaSBarry 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); 287c272547aSJed Brown #endif 288a297a907SKarl Rupp 289d25d7f95SJed Brown next->set = PETSC_TRUE; 290d25d7f95SJed Brown *((PetscInt*)next->data) = (PetscInt)lid; 291aee2cecaSBarry Smith } 292aee2cecaSBarry Smith break; 293aee2cecaSBarry Smith case OPTION_REAL: 2949566063dSJacob 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)); 2959566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 2963fc1eb6aSBarry Smith if (str[0]) { 297ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 298a4404d99SBarry Smith sscanf(str,"%e",&ir); 299570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 300570b7f6dSBarry Smith float irtemp; 301570b7f6dSBarry Smith sscanf(str,"%e",&irtemp); 302570b7f6dSBarry Smith ir = irtemp; 303ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE) 304aee2cecaSBarry Smith sscanf(str,"%le",&ir); 305ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 306d9822059SBarry Smith ir = strtoflt128(str,0); 307d9822059SBarry Smith #else 308513dbe71SLisandro Dalcin SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type"); 309a4404d99SBarry Smith #endif 310aee2cecaSBarry Smith next->set = PETSC_TRUE; 311aee2cecaSBarry Smith *((PetscReal*)next->data) = ir; 312aee2cecaSBarry Smith } 313aee2cecaSBarry Smith break; 3147781c08eSBarry Smith case OPTION_BOOL: 3159566063dSJacob 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)); 3169566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 3177781c08eSBarry Smith if (str[0]) { 3189566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringToBool(str,&bid)); 3197781c08eSBarry Smith next->set = PETSC_TRUE; 3207781c08eSBarry Smith *((PetscBool*)next->data) = bid; 3217781c08eSBarry Smith } 3227781c08eSBarry Smith break; 323aee2cecaSBarry Smith case OPTION_STRING: 3249566063dSJacob 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)); 3259566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 3263fc1eb6aSBarry Smith if (str[0]) { 327aee2cecaSBarry Smith next->set = PETSC_TRUE; 32864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3299566063dSJacob Faibussowitsch PetscCall(PetscStrdup(str,(char**)&next->data)); 3306356e834SBarry Smith } 3316356e834SBarry Smith break; 332a264d7a6SBarry Smith case OPTION_FLIST: 3339566063dSJacob Faibussowitsch PetscCall(PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data,(char*)next->data)); 3349566063dSJacob Faibussowitsch PetscCall(PetscScanString(PETSC_COMM_WORLD,512,str)); 3353cc1e11dSBarry Smith if (str[0]) { 336e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_TRUE; 3373cc1e11dSBarry Smith next->set = PETSC_TRUE; 33864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3399566063dSJacob Faibussowitsch PetscCall(PetscStrdup(str,(char**)&next->data)); 3403cc1e11dSBarry Smith } 3413cc1e11dSBarry Smith break; 342b432afdaSMatthew Knepley default: 343b432afdaSMatthew Knepley break; 3446356e834SBarry Smith } 3456356e834SBarry Smith next = next->next; 3466356e834SBarry Smith } 3476356e834SBarry Smith PetscFunctionReturn(0); 3486356e834SBarry Smith } 3496356e834SBarry Smith 350e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 351e04113cfSBarry Smith #include <petscviewersaws.h> 352d5649816SBarry Smith 353d5649816SBarry Smith static int count = 0; 354d5649816SBarry Smith 355e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void) 356d5649816SBarry Smith { 3572657e9d9SBarry Smith PetscFunctionBegin; 358d5649816SBarry Smith PetscFunctionReturn(0); 359d5649816SBarry Smith } 360d5649816SBarry Smith 3619c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n" 362a8d69d7bSBarry Smith "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n" 363a8d69d7bSBarry Smith "<script type=\"text/javascript\" src=\"https://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n" 364d1fc0251SBarry Smith "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n" 36564bbc9efSBarry Smith "<script>\n" 36664bbc9efSBarry Smith "jQuery(document).ready(function() {\n" 36764bbc9efSBarry Smith "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n" 36864bbc9efSBarry Smith "})\n" 36964bbc9efSBarry Smith "</script>\n" 37064bbc9efSBarry Smith "</head>\n"; 3711423471aSBarry Smith 3721423471aSBarry Smith /* Determines the size and style of the scroll region where PETSc options selectable from users are displayed */ 3731423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>"; 37464bbc9efSBarry Smith 375b3506946SBarry Smith /* 3767781c08eSBarry Smith PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs 377b3506946SBarry Smith 378b3506946SBarry Smith Bugs: 379b3506946SBarry Smith + All processes must traverse through the exact same set of option queries do to the call to PetscScanString() 380b3506946SBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 381b3506946SBarry Smith - Only works for PetscInt == int, PetscReal == double etc 382b3506946SBarry Smith 383b3506946SBarry Smith */ 3844416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject) 385b3506946SBarry Smith { 3864416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 387d5649816SBarry Smith static int mancount = 0; 388b3506946SBarry Smith char options[16]; 3897aab2a10SBarry Smith PetscBool changedmethod = PETSC_FALSE; 390a23eb982SSurtai Han PetscBool stopasking = PETSC_FALSE; 39188a9752cSBarry Smith char manname[16],textname[16]; 3922657e9d9SBarry Smith char dir[1024]; 393b3506946SBarry Smith 3942a409bb0SBarry Smith PetscFunctionBegin; 395b3506946SBarry Smith /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */ 396b3506946SBarry Smith sprintf(options,"Options_%d",count++); 397a297a907SKarl Rupp 3987325c4a4SBarry Smith PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */ 3991bc75a8dSBarry Smith 4009566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title")); 40183355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING)); 4029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix")); 40383355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING)); 4042657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN)); 405a23eb982SSurtai Han PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN)); 406b3506946SBarry Smith 407b3506946SBarry Smith while (next) { 408d50831c4SBarry Smith sprintf(manname,"_man_%d",mancount); 4099566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname)); 4107781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING)); 411d50831c4SBarry Smith sprintf(textname,"_text_%d",mancount++); 4129566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname)); 4137781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING)); 4149f32e415SBarry Smith 415b3506946SBarry Smith switch (next->type) { 416b3506946SBarry Smith case OPTION_HEAD: 417b3506946SBarry Smith break; 418b3506946SBarry Smith case OPTION_INT_ARRAY: 4199566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4202657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT)); 421b3506946SBarry Smith break; 422b3506946SBarry Smith case OPTION_REAL_ARRAY: 4239566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4242657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE)); 425b3506946SBarry Smith break; 426b3506946SBarry Smith case OPTION_INT: 4279566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4282657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT)); 429b3506946SBarry Smith break; 430b3506946SBarry Smith case OPTION_REAL: 4319566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4322657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE)); 433b3506946SBarry Smith break; 4347781c08eSBarry Smith case OPTION_BOOL: 4359566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4362657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN)); 4371ae3d29cSBarry Smith break; 4387781c08eSBarry Smith case OPTION_BOOL_ARRAY: 4399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4402657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN)); 44171f08665SBarry Smith break; 442b3506946SBarry Smith case OPTION_STRING: 4439566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4447781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4451ae3d29cSBarry Smith break; 4461ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 4479566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4482657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING)); 449b3506946SBarry Smith break; 450a264d7a6SBarry Smith case OPTION_FLIST: 451a264d7a6SBarry Smith { 452a264d7a6SBarry Smith PetscInt ntext; 4539566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4547781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4559566063dSJacob Faibussowitsch PetscCall(PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext)); 456a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 457a264d7a6SBarry Smith } 458a264d7a6SBarry Smith break; 4591ae3d29cSBarry Smith case OPTION_ELIST: 460a264d7a6SBarry Smith { 461a264d7a6SBarry Smith PetscInt ntext = next->nlist; 4629566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 4637781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((ntext+1),(char***)&next->edata)); 4659566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(next->edata,next->list,ntext*sizeof(char*))); 466a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 467a264d7a6SBarry Smith } 468a264d7a6SBarry Smith break; 469b3506946SBarry Smith default: 470b3506946SBarry Smith break; 471b3506946SBarry Smith } 472b3506946SBarry Smith next = next->next; 473b3506946SBarry Smith } 474b3506946SBarry Smith 475b3506946SBarry Smith /* wait until accessor has unlocked the memory */ 47664bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader)); 47764bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom)); 4789566063dSJacob Faibussowitsch PetscCall(PetscSAWsBlock()); 47964bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Header,("index.html")); 48064bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2)); 481b3506946SBarry Smith 48288a9752cSBarry Smith /* determine if any values have been set in GUI */ 48383355fc5SBarry Smith next = PetscOptionsObject->next; 48488a9752cSBarry Smith while (next) { 4859566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option)); 486f7b25cf6SBarry Smith PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set)); 48788a9752cSBarry Smith next = next->next; 48888a9752cSBarry Smith } 48988a9752cSBarry Smith 490b3506946SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 491f7b25cf6SBarry Smith if (changedmethod) PetscOptionsObject->count = -2; 492b3506946SBarry Smith 493a23eb982SSurtai Han if (stopasking) { 494a23eb982SSurtai Han PetscOptionsPublish = PETSC_FALSE; 49512655325SBarry Smith PetscOptionsObject->count = 0;//do not ask for same thing again 496a23eb982SSurtai Han } 497a23eb982SSurtai Han 4989a492a5cSBarry Smith PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options")); 499b3506946SBarry Smith PetscFunctionReturn(0); 500b3506946SBarry Smith } 501b3506946SBarry Smith #endif 502b3506946SBarry Smith 5034416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject) 50453acd3b1SBarry Smith { 5054416b707SBarry Smith PetscOptionItem last; 5066356e834SBarry Smith char option[256],value[1024],tmp[32]; 507330cf3c9SBarry Smith size_t j; 50853acd3b1SBarry Smith 50953acd3b1SBarry Smith PetscFunctionBegin; 51083355fc5SBarry Smith if (PetscOptionsObject->next) { 51183355fc5SBarry Smith if (!PetscOptionsObject->count) { 512a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS) 5139566063dSJacob Faibussowitsch PetscCall(PetscOptionsSAWsInput(PetscOptionsObject)); 514b3506946SBarry Smith #else 5159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetFromTextInput(PetscOptionsObject)); 516b3506946SBarry Smith #endif 517aee2cecaSBarry Smith } 518aee2cecaSBarry Smith } 5196356e834SBarry Smith 5209566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->title)); 5216356e834SBarry Smith 522e26ddf31SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 523e55864a3SBarry Smith if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2; 5247a72a596SBarry Smith /* reset alreadyprinted flag */ 525e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = PETSC_FALSE; 526e55864a3SBarry Smith if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE; 527e55864a3SBarry Smith PetscOptionsObject->object = NULL; 52853acd3b1SBarry Smith 529e55864a3SBarry Smith while (PetscOptionsObject->next) { 530e55864a3SBarry Smith if (PetscOptionsObject->next->set) { 531e55864a3SBarry Smith if (PetscOptionsObject->prefix) { 5329566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(option,"-")); 5339566063dSJacob Faibussowitsch PetscCall(PetscStrcat(option,PetscOptionsObject->prefix)); 5349566063dSJacob Faibussowitsch PetscCall(PetscStrcat(option,PetscOptionsObject->next->option+1)); 5359566063dSJacob Faibussowitsch } else PetscCall(PetscStrcpy(option,PetscOptionsObject->next->option)); 5366356e834SBarry Smith 537e55864a3SBarry Smith switch (PetscOptionsObject->next->type) { 5386356e834SBarry Smith case OPTION_HEAD: 5396356e834SBarry Smith break; 5406356e834SBarry Smith case OPTION_INT_ARRAY: 541e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]); 542e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 543e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]); 5449566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,",")); 5459566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,tmp)); 5466356e834SBarry Smith } 5476356e834SBarry Smith break; 5486356e834SBarry Smith case OPTION_INT: 549e55864a3SBarry Smith sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data); 5506356e834SBarry Smith break; 5516356e834SBarry Smith case OPTION_REAL: 552e55864a3SBarry Smith sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data); 5536356e834SBarry Smith break; 5546356e834SBarry Smith case OPTION_REAL_ARRAY: 555e55864a3SBarry Smith sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]); 556e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 557e55864a3SBarry Smith sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]); 5589566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,",")); 5599566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,tmp)); 5606356e834SBarry Smith } 5616356e834SBarry Smith break; 562050cccc3SHong Zhang case OPTION_SCALAR_ARRAY: 56395f3a755SJose E. Roman sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0])); 564050cccc3SHong Zhang for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 56595f3a755SJose E. Roman sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j])); 5669566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,",")); 5679566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,tmp)); 568050cccc3SHong Zhang } 569050cccc3SHong Zhang break; 5707781c08eSBarry Smith case OPTION_BOOL: 571e55864a3SBarry Smith sprintf(value,"%d",*(int*)PetscOptionsObject->next->data); 5726356e834SBarry Smith break; 5737781c08eSBarry Smith case OPTION_BOOL_ARRAY: 574e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]); 575e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 576e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]); 5779566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,",")); 5789566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,tmp)); 5791ae3d29cSBarry Smith } 5801ae3d29cSBarry Smith break; 581a264d7a6SBarry Smith case OPTION_FLIST: 5829566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(value,(char*)PetscOptionsObject->next->data)); 5836991f827SBarry Smith break; 5841ae3d29cSBarry Smith case OPTION_ELIST: 5859566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(value,(char*)PetscOptionsObject->next->data)); 5866356e834SBarry Smith break; 5871ae3d29cSBarry Smith case OPTION_STRING: 5889566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(value,(char*)PetscOptionsObject->next->data)); 589d51da6bfSBarry Smith break; 5901ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 591e55864a3SBarry Smith sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]); 592e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 593e55864a3SBarry Smith sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]); 5949566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,",")); 5959566063dSJacob Faibussowitsch PetscCall(PetscStrcat(value,tmp)); 5961ae3d29cSBarry Smith } 5976356e834SBarry Smith break; 5986356e834SBarry Smith } 5999566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(PetscOptionsObject->options,option,value)); 6006356e834SBarry Smith } 6016991f827SBarry Smith if (PetscOptionsObject->next->type == OPTION_ELIST) { 6029566063dSJacob Faibussowitsch PetscCall(PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list)); 6036991f827SBarry Smith } 6049566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->text)); 6059566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->option)); 6069566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->man)); 6079566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->edata)); 608c979a496SBarry Smith 60983355fc5SBarry Smith if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)) { 61083355fc5SBarry Smith free(PetscOptionsObject->next->data); 611c979a496SBarry Smith } else { 6129566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->next->data)); 613c979a496SBarry Smith } 6147781c08eSBarry Smith 61583355fc5SBarry Smith last = PetscOptionsObject->next; 61683355fc5SBarry Smith PetscOptionsObject->next = PetscOptionsObject->next->next; 6179566063dSJacob Faibussowitsch PetscCall(PetscFree(last)); 6186356e834SBarry Smith } 6199566063dSJacob Faibussowitsch PetscCall(PetscFree(PetscOptionsObject->prefix)); 62002c9f0b5SLisandro Dalcin PetscOptionsObject->next = NULL; 62153acd3b1SBarry Smith PetscFunctionReturn(0); 62253acd3b1SBarry Smith } 62353acd3b1SBarry Smith 62488aa4217SBarry Smith /*MC 62553acd3b1SBarry Smith PetscOptionsEnum - Gets the enum value for a particular option in the database. 62653acd3b1SBarry Smith 6273f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 62853acd3b1SBarry Smith 62988aa4217SBarry Smith Synopsis: 63088aa4217SBarry Smith #include "petscsys.h" 6313a89f35bSSatish Balay PetscErrorCode PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool *set) 63288aa4217SBarry Smith 63353acd3b1SBarry Smith Input Parameters: 63453acd3b1SBarry Smith + opt - option name 63553acd3b1SBarry Smith . text - short string that describes the option 63653acd3b1SBarry Smith . man - manual page with additional information on option 63753acd3b1SBarry Smith . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 6380fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6390fdccdaeSBarry Smith $ PetscOptionsEnum(..., obj->value,&object->value,...) or 6400fdccdaeSBarry Smith $ value = defaultvalue 6410fdccdaeSBarry Smith $ PetscOptionsEnum(..., value,&value,&flg); 6420fdccdaeSBarry Smith $ if (flg) { 64353acd3b1SBarry Smith 644d8d19677SJose E. Roman Output Parameters: 64553acd3b1SBarry Smith + value - the value to return 646b32e0204SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 64753acd3b1SBarry Smith 64853acd3b1SBarry Smith Level: beginner 64953acd3b1SBarry Smith 65095452b02SPatrick Sanan Notes: 65195452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 65253acd3b1SBarry Smith 65353acd3b1SBarry Smith list is usually something like PCASMTypes or some other predefined list of enum names 65453acd3b1SBarry Smith 6552efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 6562efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 6572efd9cb1SBarry Smith 658989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 659989712b9SBarry Smith 660db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 661db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 662db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 663db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 664*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 665db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 666db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 66788aa4217SBarry Smith M*/ 66888aa4217SBarry Smith 6694416b707SBarry Smith PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool *set) 67053acd3b1SBarry Smith { 67153acd3b1SBarry Smith PetscInt ntext = 0; 672aa5bb8c0SSatish Balay PetscInt tval; 673ace3abfcSBarry Smith PetscBool tflg; 67453acd3b1SBarry Smith 67553acd3b1SBarry Smith PetscFunctionBegin; 67653acd3b1SBarry Smith while (list[ntext++]) { 67708401ef6SPierre Jolivet PetscCheck(ntext <= 50,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 67853acd3b1SBarry Smith } 67908401ef6SPierre Jolivet PetscCheck(ntext >= 3,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 68053acd3b1SBarry Smith ntext -= 3; 6819566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg)); 682aa5bb8c0SSatish Balay /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 683aa5bb8c0SSatish Balay if (tflg) *value = (PetscEnum)tval; 684aa5bb8c0SSatish Balay if (set) *set = tflg; 68553acd3b1SBarry Smith PetscFunctionReturn(0); 68653acd3b1SBarry Smith } 68753acd3b1SBarry Smith 68888aa4217SBarry Smith /*MC 689d3e47460SLisandro Dalcin PetscOptionsEnumArray - Gets an array of enum values for a particular 690d3e47460SLisandro Dalcin option in the database. 691d3e47460SLisandro Dalcin 692d3e47460SLisandro Dalcin Logically Collective on the communicator passed in PetscOptionsBegin() 693d3e47460SLisandro Dalcin 69488aa4217SBarry Smith Synopsis: 69588aa4217SBarry Smith #include "petscsys.h" 6963a89f35bSSatish Balay PetscErrorCode PetscOptionsEnumArray(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool *set) 69788aa4217SBarry Smith 698d3e47460SLisandro Dalcin Input Parameters: 699d3e47460SLisandro Dalcin + opt - the option one is seeking 700d3e47460SLisandro Dalcin . text - short string describing option 701d3e47460SLisandro Dalcin . man - manual page for option 70222d58ec6SMatthew G. Knepley . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 703c89788bdSBarry Smith - n - maximum number of values allowed in the value array 704d3e47460SLisandro Dalcin 705d8d19677SJose E. Roman Output Parameters: 706d3e47460SLisandro Dalcin + value - location to copy values 707d3e47460SLisandro Dalcin . n - actual number of values found 708d3e47460SLisandro Dalcin - set - PETSC_TRUE if found, else PETSC_FALSE 709d3e47460SLisandro Dalcin 710d3e47460SLisandro Dalcin Level: beginner 711d3e47460SLisandro Dalcin 712d3e47460SLisandro Dalcin Notes: 713d3e47460SLisandro Dalcin The array must be passed as a comma separated list. 714d3e47460SLisandro Dalcin 715d3e47460SLisandro Dalcin There must be no intervening spaces between the values. 716d3e47460SLisandro Dalcin 717d3e47460SLisandro Dalcin Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 718d3e47460SLisandro Dalcin 719db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 720db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 721db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 722*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 723db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 724db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRealArray()` 72588aa4217SBarry Smith M*/ 72688aa4217SBarry Smith 7274416b707SBarry Smith PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool *set) 728d3e47460SLisandro Dalcin { 729d3e47460SLisandro Dalcin PetscInt i,nlist = 0; 7304416b707SBarry Smith PetscOptionItem amsopt; 731d3e47460SLisandro Dalcin 732d3e47460SLisandro Dalcin PetscFunctionBegin; 73308401ef6SPierre 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"); 73408401ef6SPierre Jolivet PetscCheck(nlist >= 3,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 735d3e47460SLisandro Dalcin nlist -= 3; /* drop enum name, prefix, and null termination */ 736d3e47460SLisandro Dalcin if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */ 737d3e47460SLisandro Dalcin PetscEnum *vals; 7389566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt)); 7399566063dSJacob Faibussowitsch PetscCall(PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list)); 740d3e47460SLisandro Dalcin amsopt->nlist = nlist; 7419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n,(PetscEnum**)&amsopt->data)); 742d3e47460SLisandro Dalcin amsopt->arraylength = *n; 743d3e47460SLisandro Dalcin vals = (PetscEnum*)amsopt->data; 744d3e47460SLisandro Dalcin for (i=0; i<*n; i++) vals[i] = value[i]; 745d3e47460SLisandro Dalcin } 7469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set)); 747d3e47460SLisandro Dalcin if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 7489566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]])); 7499566063dSJacob Faibussowitsch for (i=1; i<*n; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]])); 7509566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text)); 7519566063dSJacob Faibussowitsch for (i=0; i<nlist; i++) PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i])); 7529566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man))); 753d3e47460SLisandro Dalcin } 754d3e47460SLisandro Dalcin PetscFunctionReturn(0); 755d3e47460SLisandro Dalcin } 756d3e47460SLisandro Dalcin 75788aa4217SBarry Smith /*MC 7585a856986SBarry Smith PetscOptionsBoundedInt - Gets an integer value greater than or equal a given bound for a particular option in the database. 7595a856986SBarry Smith 7605a856986SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 7615a856986SBarry Smith 76288aa4217SBarry Smith Synopsis: 76388aa4217SBarry Smith #include "petscsys.h" 7643a89f35bSSatish Balay PetscErrorCode PetscOptionsBoundInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg,PetscInt bound) 76588aa4217SBarry Smith 7665a856986SBarry Smith Input Parameters: 7675a856986SBarry Smith + opt - option name 7685a856986SBarry Smith . text - short string that describes the option 7695a856986SBarry Smith . man - manual page with additional information on option 7705a856986SBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 77188aa4217SBarry Smith $ PetscOptionsInt(..., obj->value,&obj->value,...) or 7725a856986SBarry Smith $ value = defaultvalue 7735a856986SBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 7745a856986SBarry Smith $ if (flg) { 77588aa4217SBarry Smith - bound - the requested value should be greater than or equal this bound or an error is generated 7765a856986SBarry Smith 777d8d19677SJose E. Roman Output Parameters: 7785a856986SBarry Smith + value - the integer value to return 7795a856986SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 7805a856986SBarry Smith 7815a856986SBarry Smith Notes: 7825a856986SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 7835a856986SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 7845a856986SBarry Smith 7855a856986SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 7865a856986SBarry Smith 7875a856986SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 7885a856986SBarry Smith 7895a856986SBarry Smith Level: beginner 7905a856986SBarry Smith 791db781477SPatrick Sanan .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 792db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 793db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 794db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 795*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 796db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 797db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 79888aa4217SBarry Smith M*/ 7995a856986SBarry Smith 80088aa4217SBarry Smith /*MC 8015a856986SBarry Smith PetscOptionsRangeInt - Gets an integer value within a range of values for a particular option in the database. 8025a856986SBarry Smith 8035a856986SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 8045a856986SBarry Smith 80588aa4217SBarry Smith Synopsis: 80688aa4217SBarry Smith #include "petscsys.h" 8073a89f35bSSatish Balay PetscErrorCode PetscOptionsRangeInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg,PetscInt lb,PetscInt ub) 80888aa4217SBarry Smith 8095a856986SBarry Smith Input Parameters: 8105a856986SBarry Smith + opt - option name 8115a856986SBarry Smith . text - short string that describes the option 8125a856986SBarry Smith . man - manual page with additional information on option 8135a856986SBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 81488aa4217SBarry Smith $ PetscOptionsInt(..., obj->value,&obj->value,...) or 8155a856986SBarry Smith $ value = defaultvalue 8165a856986SBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 8175a856986SBarry Smith $ if (flg) { 81888aa4217SBarry Smith . lb - the lower bound, provided value must be greater than or equal to this value or an error is generated 81988aa4217SBarry Smith - ub - the upper bound, provided value must be less than or equal to this value or an error is generated 8205a856986SBarry Smith 821d8d19677SJose E. Roman Output Parameters: 8225a856986SBarry Smith + value - the integer value to return 8235a856986SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 8245a856986SBarry Smith 8255a856986SBarry Smith Notes: 8265a856986SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 8275a856986SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 8285a856986SBarry Smith 8295a856986SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 8305a856986SBarry Smith 8315a856986SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 8325a856986SBarry Smith 8335a856986SBarry Smith Level: beginner 8345a856986SBarry Smith 835db781477SPatrick Sanan .seealso: `PetscOptionsInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 836db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsBoundedInt()` 837db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 838db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 839*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 840db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 841db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 84288aa4217SBarry Smith M*/ 8435a856986SBarry Smith 84488aa4217SBarry Smith /*MC 84553acd3b1SBarry Smith PetscOptionsInt - Gets the integer value for a particular option in the database. 84653acd3b1SBarry Smith 8473f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 84853acd3b1SBarry Smith 84988aa4217SBarry Smith Synopsis: 85088aa4217SBarry Smith #include "petscsys.h" 8516eeb591dSVaclav Hapla PetscErrorCode PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *flg)) 85288aa4217SBarry Smith 85353acd3b1SBarry Smith Input Parameters: 85453acd3b1SBarry Smith + opt - option name 85553acd3b1SBarry Smith . text - short string that describes the option 85653acd3b1SBarry Smith . man - manual page with additional information on option 8570fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 85888aa4217SBarry Smith $ PetscOptionsInt(..., obj->value,&obj->value,...) or 8590fdccdaeSBarry Smith $ value = defaultvalue 8600fdccdaeSBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 8610fdccdaeSBarry Smith $ if (flg) { 86253acd3b1SBarry Smith 863d8d19677SJose E. Roman Output Parameters: 86453acd3b1SBarry Smith + value - the integer value to return 86553acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 86653acd3b1SBarry Smith 86795452b02SPatrick Sanan Notes: 86895452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 8692efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 8702efd9cb1SBarry Smith 871989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 872989712b9SBarry Smith 87395452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 87453acd3b1SBarry Smith 8755a856986SBarry Smith Level: beginner 8765a856986SBarry Smith 877db781477SPatrick Sanan .seealso: `PetscOptionsBoundedInt()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 878db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, `PetscOptionsRangeInt()` 879db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 880db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 881*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 882db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 883db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 88488aa4217SBarry Smith M*/ 88588aa4217SBarry Smith 8865a856986SBarry Smith PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *set,PetscInt lb,PetscInt ub) 88753acd3b1SBarry Smith { 8884416b707SBarry Smith PetscOptionItem amsopt; 88912655325SBarry Smith PetscBool wasset; 89053acd3b1SBarry Smith 89153acd3b1SBarry Smith PetscFunctionBegin; 89208401ef6SPierre Jolivet PetscCheck(currentvalue >= lb,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Current value %" PetscInt_FMT " less than allowed bound %" PetscInt_FMT,currentvalue,lb); 89308401ef6SPierre Jolivet PetscCheck(currentvalue <= ub,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Current value %" PetscInt_FMT " greater than allowed bound %" PetscInt_FMT,currentvalue,ub); 894e55864a3SBarry Smith if (!PetscOptionsObject->count) { 8959566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt)); 8969566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscInt),&amsopt->data)); 89712655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 8983e211508SBarry Smith 8999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,¤tvalue,&wasset)); 9003e211508SBarry Smith if (wasset) { 90112655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 9023e211508SBarry Smith } 903af6d86caSBarry Smith } 9049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&wasset)); 905cc73adaaSBarry 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); 906cc73adaaSBarry 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); 90744ef3d73SBarry Smith if (set) *set = wasset; 908e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9099566063dSJacob 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))); 91053acd3b1SBarry Smith } 91153acd3b1SBarry Smith PetscFunctionReturn(0); 91253acd3b1SBarry Smith } 91353acd3b1SBarry Smith 91488aa4217SBarry Smith /*MC 91553acd3b1SBarry Smith PetscOptionsString - Gets the string value for a particular option in the database. 91653acd3b1SBarry Smith 9173f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 91853acd3b1SBarry Smith 91988aa4217SBarry Smith Synopsis: 92088aa4217SBarry Smith #include "petscsys.h" 9213a89f35bSSatish Balay PetscErrorCode PetscOptionsString(const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool *set) 92288aa4217SBarry Smith 92353acd3b1SBarry Smith Input Parameters: 92453acd3b1SBarry Smith + opt - option name 92553acd3b1SBarry Smith . text - short string that describes the option 92653acd3b1SBarry Smith . man - manual page with additional information on option 9270fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 928bcbf2dc5SJed Brown - len - length of the result string including null terminator 92953acd3b1SBarry Smith 930d8d19677SJose E. Roman Output Parameters: 93153acd3b1SBarry Smith + value - the value to return 93253acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 93353acd3b1SBarry Smith 93453acd3b1SBarry Smith Level: beginner 93553acd3b1SBarry Smith 93695452b02SPatrick Sanan Notes: 93795452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 93853acd3b1SBarry Smith 9397fccdfe4SBarry 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). 9407fccdfe4SBarry Smith 9412efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 9422efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 9432efd9cb1SBarry Smith 944989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 945989712b9SBarry Smith 946db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 947db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 948db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 949db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 950*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 951db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 952db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 95388aa4217SBarry Smith M*/ 95488aa4217SBarry Smith 9554416b707SBarry Smith PetscErrorCode PetscOptionsString_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool *set) 95653acd3b1SBarry Smith { 9574416b707SBarry Smith PetscOptionItem amsopt; 95844ef3d73SBarry Smith PetscBool lset; 95953acd3b1SBarry Smith 96053acd3b1SBarry Smith PetscFunctionBegin; 9611a1499c8SBarry Smith if (!PetscOptionsObject->count) { 9629566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt)); 96364facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 9649566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data)); 965af6d86caSBarry Smith } 9669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset)); 96744ef3d73SBarry Smith if (set) *set = lset; 968e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9699566063dSJacob 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))); 97053acd3b1SBarry Smith } 97153acd3b1SBarry Smith PetscFunctionReturn(0); 97253acd3b1SBarry Smith } 97353acd3b1SBarry Smith 97488aa4217SBarry Smith /*MC 97553acd3b1SBarry Smith PetscOptionsReal - Gets the PetscReal value for a particular option in the database. 97653acd3b1SBarry Smith 9773f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 97853acd3b1SBarry Smith 97988aa4217SBarry Smith Synopsis: 98088aa4217SBarry Smith #include "petscsys.h" 9813a89f35bSSatish Balay PetscErrorCode PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 98288aa4217SBarry Smith 98353acd3b1SBarry Smith Input Parameters: 98453acd3b1SBarry Smith + opt - option name 98553acd3b1SBarry Smith . text - short string that describes the option 98653acd3b1SBarry Smith . man - manual page with additional information on option 9870fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 98888aa4217SBarry Smith $ PetscOptionsReal(..., obj->value,&obj->value,...) or 9890fdccdaeSBarry Smith $ value = defaultvalue 9900fdccdaeSBarry Smith $ PetscOptionsReal(..., value,&value,&flg); 9910fdccdaeSBarry Smith $ if (flg) { 99253acd3b1SBarry Smith 993d8d19677SJose E. Roman Output Parameters: 99453acd3b1SBarry Smith + value - the value to return 99553acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 99653acd3b1SBarry Smith 99795452b02SPatrick Sanan Notes: 99895452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 9992efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 10002efd9cb1SBarry Smith 1001989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 1002989712b9SBarry Smith 100395452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 100453acd3b1SBarry Smith 10055a856986SBarry Smith Level: beginner 10065a856986SBarry Smith 1007db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1008db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1009db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1010db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1011*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1012db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1013db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 101488aa4217SBarry Smith M*/ 101588aa4217SBarry Smith 10164416b707SBarry Smith PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 101753acd3b1SBarry Smith { 10184416b707SBarry Smith PetscOptionItem amsopt; 101944ef3d73SBarry Smith PetscBool lset; 102053acd3b1SBarry Smith 102153acd3b1SBarry Smith PetscFunctionBegin; 1022e55864a3SBarry Smith if (!PetscOptionsObject->count) { 10239566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt)); 10249566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscReal),&amsopt->data)); 1025a297a907SKarl Rupp 10260fdccdaeSBarry Smith *(PetscReal*)amsopt->data = currentvalue; 1027538aa990SBarry Smith } 10289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,&lset)); 102944ef3d73SBarry Smith if (set) *set = lset; 10301a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10319566063dSJacob 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))); 103253acd3b1SBarry Smith } 103353acd3b1SBarry Smith PetscFunctionReturn(0); 103453acd3b1SBarry Smith } 103553acd3b1SBarry Smith 103688aa4217SBarry Smith /*MC 103753acd3b1SBarry Smith PetscOptionsScalar - Gets the scalar value for a particular option in the database. 103853acd3b1SBarry Smith 10393f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 104053acd3b1SBarry Smith 104188aa4217SBarry Smith Synopsis: 104288aa4217SBarry Smith #include "petscsys.h" 10433a89f35bSSatish Balay PetscErrorCode PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 104488aa4217SBarry Smith 104553acd3b1SBarry Smith Input Parameters: 104653acd3b1SBarry Smith + opt - option name 104753acd3b1SBarry Smith . text - short string that describes the option 104853acd3b1SBarry Smith . man - manual page with additional information on option 10490fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 105088aa4217SBarry Smith $ PetscOptionsScalar(..., obj->value,&obj->value,...) or 10510fdccdaeSBarry Smith $ value = defaultvalue 10520fdccdaeSBarry Smith $ PetscOptionsScalar(..., value,&value,&flg); 10530fdccdaeSBarry Smith $ if (flg) { 10540fdccdaeSBarry Smith 1055d8d19677SJose E. Roman Output Parameters: 105653acd3b1SBarry Smith + value - the value to return 105753acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 105853acd3b1SBarry Smith 105995452b02SPatrick Sanan Notes: 106095452b02SPatrick Sanan If the user does not supply the option at all value is NOT changed. Thus 10612efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 10622efd9cb1SBarry Smith 1063989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 1064989712b9SBarry Smith 106595452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 106653acd3b1SBarry Smith 10675a856986SBarry Smith Level: beginner 10685a856986SBarry Smith 1069db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1070db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1071db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1072db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1073*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1074db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1075db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 107688aa4217SBarry Smith M*/ 107788aa4217SBarry Smith 10784416b707SBarry Smith PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 107953acd3b1SBarry Smith { 108053acd3b1SBarry Smith PetscFunctionBegin; 108153acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 10829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal(opt,text,man,currentvalue,value,set)); 108353acd3b1SBarry Smith #else 10849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set)); 108553acd3b1SBarry Smith #endif 108653acd3b1SBarry Smith PetscFunctionReturn(0); 108753acd3b1SBarry Smith } 108853acd3b1SBarry Smith 108988aa4217SBarry Smith /*MC 109090d69ab7SBarry 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 109190d69ab7SBarry Smith its value is set to false. 109253acd3b1SBarry Smith 10933f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 109453acd3b1SBarry Smith 109588aa4217SBarry Smith Synopsis: 109688aa4217SBarry Smith #include "petscsys.h" 10973a89f35bSSatish Balay PetscErrorCode PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool *flg) 109888aa4217SBarry Smith 109953acd3b1SBarry Smith Input Parameters: 110053acd3b1SBarry Smith + opt - option name 110153acd3b1SBarry Smith . text - short string that describes the option 110253acd3b1SBarry Smith - man - manual page with additional information on option 110353acd3b1SBarry Smith 110453acd3b1SBarry Smith Output Parameter: 110553acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 110653acd3b1SBarry Smith 110753acd3b1SBarry Smith Level: beginner 110853acd3b1SBarry Smith 110995452b02SPatrick Sanan Notes: 111095452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 111153acd3b1SBarry Smith 1112db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1113db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1114db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1115db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1116*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1117db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1118db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 111988aa4217SBarry Smith M*/ 112088aa4217SBarry Smith 11214416b707SBarry Smith PetscErrorCode PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 112253acd3b1SBarry Smith { 11234416b707SBarry Smith PetscOptionItem amsopt; 112453acd3b1SBarry Smith 112553acd3b1SBarry Smith PetscFunctionBegin; 1126e55864a3SBarry Smith if (!PetscOptionsObject->count) { 11279566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt)); 11289566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool),&amsopt->data)); 1129a297a907SKarl Rupp 1130ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 11311ae3d29cSBarry Smith } 11329566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg)); 1133e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 11349566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man))); 113553acd3b1SBarry Smith } 113653acd3b1SBarry Smith PetscFunctionReturn(0); 113753acd3b1SBarry Smith } 113853acd3b1SBarry Smith 113988aa4217SBarry Smith /*MC 1140a264d7a6SBarry Smith PetscOptionsFList - Puts a list of option values that a single one may be selected from 114153acd3b1SBarry Smith 11423f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 114353acd3b1SBarry Smith 114488aa4217SBarry Smith Synopsis: 114588aa4217SBarry Smith #include "petscsys.h" 11463a89f35bSSatish Balay PetscErrorCode PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool *set) 114788aa4217SBarry Smith 114853acd3b1SBarry Smith Input Parameters: 114953acd3b1SBarry Smith + opt - option name 115053acd3b1SBarry Smith . text - short string that describes the option 115153acd3b1SBarry Smith . man - manual page with additional information on option 115253acd3b1SBarry Smith . list - the possible choices 11530fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 11540fdccdaeSBarry Smith $ PetscOptionsFlist(..., obj->value,value,len,&flg); 11550fdccdaeSBarry Smith $ if (flg) { 11563cc1e11dSBarry Smith - len - the length of the character array value 115753acd3b1SBarry Smith 1158d8d19677SJose E. Roman Output Parameters: 115953acd3b1SBarry Smith + value - the value to return 116053acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 116153acd3b1SBarry Smith 116253acd3b1SBarry Smith Level: intermediate 116353acd3b1SBarry Smith 116495452b02SPatrick Sanan Notes: 116595452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 116653acd3b1SBarry Smith 11672efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 11682efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 11692efd9cb1SBarry Smith 1170989712b9SBarry Smith The default/currentvalue passed into this routine does not get transferred to the output value variable automatically. 1171989712b9SBarry Smith 117253acd3b1SBarry Smith See PetscOptionsEList() for when the choices are given in a string array 117353acd3b1SBarry Smith 117453acd3b1SBarry Smith To get a listing of all currently specified options, 117588c29154SBarry Smith see PetscOptionsView() or PetscOptionsGetAll() 117653acd3b1SBarry Smith 1177eabe10d7SBarry Smith Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list 1178eabe10d7SBarry Smith 1179db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1180db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1181db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1182*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1183db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1184db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsEnum()` 118588aa4217SBarry Smith M*/ 118688aa4217SBarry Smith 11874416b707SBarry Smith 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) 118853acd3b1SBarry Smith { 11894416b707SBarry Smith PetscOptionItem amsopt; 119044ef3d73SBarry Smith PetscBool lset; 119153acd3b1SBarry Smith 119253acd3b1SBarry Smith PetscFunctionBegin; 11931a1499c8SBarry Smith if (!PetscOptionsObject->count) { 11949566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt)); 119564facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 11969566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data)); 11973cc1e11dSBarry Smith amsopt->flist = list; 11983cc1e11dSBarry Smith } 11999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,&lset)); 120044ef3d73SBarry Smith if (set) *set = lset; 12011a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 12029566063dSJacob Faibussowitsch PetscCall(PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue,lset && value ? value : currentvalue)); 120353acd3b1SBarry Smith } 120453acd3b1SBarry Smith PetscFunctionReturn(0); 120553acd3b1SBarry Smith } 120653acd3b1SBarry Smith 120788aa4217SBarry Smith /*MC 120853acd3b1SBarry Smith PetscOptionsEList - Puts a list of option values that a single one may be selected from 120953acd3b1SBarry Smith 12103f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 121153acd3b1SBarry Smith 121288aa4217SBarry Smith Synopsis: 121388aa4217SBarry Smith #include "petscsys.h" 12143a89f35bSSatish 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) 121588aa4217SBarry Smith 121653acd3b1SBarry Smith Input Parameters: 121753acd3b1SBarry Smith + opt - option name 121853acd3b1SBarry Smith . ltext - short string that describes the option 121953acd3b1SBarry Smith . man - manual page with additional information on option 1220a264d7a6SBarry Smith . list - the possible choices (one of these must be selected, anything else is invalid) 122153acd3b1SBarry Smith . ntext - number of choices 12220fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 12230fdccdaeSBarry Smith $ PetscOptionsElist(..., obj->value,&value,&flg); 12240fdccdaeSBarry Smith $ if (flg) { 12250fdccdaeSBarry Smith 1226d8d19677SJose E. Roman Output Parameters: 122753acd3b1SBarry Smith + value - the index of the value to return 122853acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 122953acd3b1SBarry Smith 123053acd3b1SBarry Smith Level: intermediate 123153acd3b1SBarry Smith 123295452b02SPatrick Sanan Notes: 123395452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 123453acd3b1SBarry Smith 12352efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 12362efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 12372efd9cb1SBarry Smith 1238a264d7a6SBarry Smith See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 123953acd3b1SBarry Smith 1240db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1241db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1242db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1243*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1244db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1245db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEnum()` 124688aa4217SBarry Smith M*/ 124788aa4217SBarry Smith 12484416b707SBarry Smith 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) 124953acd3b1SBarry Smith { 125053acd3b1SBarry Smith PetscInt i; 12514416b707SBarry Smith PetscOptionItem amsopt; 125244ef3d73SBarry Smith PetscBool lset; 125353acd3b1SBarry Smith 125453acd3b1SBarry Smith PetscFunctionBegin; 12551a1499c8SBarry Smith if (!PetscOptionsObject->count) { 12569566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt)); 125764facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 12589566063dSJacob Faibussowitsch PetscCall(PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data)); 12599566063dSJacob Faibussowitsch PetscCall(PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list)); 12601ae3d29cSBarry Smith amsopt->nlist = ntext; 12611ae3d29cSBarry Smith } 12629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,&lset)); 126344ef3d73SBarry Smith if (set) *set = lset; 12641a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 12659566063dSJacob 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)); 126653acd3b1SBarry Smith for (i=0; i<ntext; i++) { 12679566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i])); 126853acd3b1SBarry Smith } 12699566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man))); 127053acd3b1SBarry Smith } 127153acd3b1SBarry Smith PetscFunctionReturn(0); 127253acd3b1SBarry Smith } 127353acd3b1SBarry Smith 127488aa4217SBarry Smith /*MC 1275acfcf0e5SJed Brown PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 1276d5649816SBarry Smith which at most a single value can be true. 127753acd3b1SBarry Smith 12783f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 127953acd3b1SBarry Smith 128088aa4217SBarry Smith Synopsis: 128188aa4217SBarry Smith #include "petscsys.h" 12823a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool *flg) 128388aa4217SBarry Smith 128453acd3b1SBarry Smith Input Parameters: 128553acd3b1SBarry Smith + opt - option name 128653acd3b1SBarry Smith . text - short string that describes the option 128753acd3b1SBarry Smith - man - manual page with additional information on option 128853acd3b1SBarry Smith 128953acd3b1SBarry Smith Output Parameter: 129053acd3b1SBarry Smith . flg - whether that option was set or not 129153acd3b1SBarry Smith 129253acd3b1SBarry Smith Level: intermediate 129353acd3b1SBarry Smith 129495452b02SPatrick Sanan Notes: 129595452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 129653acd3b1SBarry Smith 1297acfcf0e5SJed Brown Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd() 129853acd3b1SBarry Smith 1299db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1300db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1301db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1302*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1303db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1304db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 130588aa4217SBarry Smith M*/ 130688aa4217SBarry Smith 13074416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 130853acd3b1SBarry Smith { 13094416b707SBarry Smith PetscOptionItem amsopt; 131053acd3b1SBarry Smith 131153acd3b1SBarry Smith PetscFunctionBegin; 1312e55864a3SBarry Smith if (!PetscOptionsObject->count) { 13139566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt)); 13149566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool),&amsopt->data)); 1315a297a907SKarl Rupp 1316ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 13171ae3d29cSBarry Smith } 131868b16fdaSBarry Smith *flg = PETSC_FALSE; 13199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL)); 1320e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 13219566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," Pick at most one of -------------\n")); 13229566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man))); 132353acd3b1SBarry Smith } 132453acd3b1SBarry Smith PetscFunctionReturn(0); 132553acd3b1SBarry Smith } 132653acd3b1SBarry Smith 132788aa4217SBarry Smith /*MC 1328acfcf0e5SJed Brown PetscOptionsBoolGroup - One in a series of logical queries on the options database for 1329d5649816SBarry Smith which at most a single value can be true. 133053acd3b1SBarry Smith 13313f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 133253acd3b1SBarry Smith 133388aa4217SBarry Smith Synopsis: 133488aa4217SBarry Smith #include "petscsys.h" 13353a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool *flg) 133688aa4217SBarry Smith 133753acd3b1SBarry Smith Input Parameters: 133853acd3b1SBarry Smith + opt - option name 133953acd3b1SBarry Smith . text - short string that describes the option 134053acd3b1SBarry Smith - man - manual page with additional information on option 134153acd3b1SBarry Smith 134253acd3b1SBarry Smith Output Parameter: 134353acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 134453acd3b1SBarry Smith 134553acd3b1SBarry Smith Level: intermediate 134653acd3b1SBarry Smith 134795452b02SPatrick Sanan Notes: 134895452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 134953acd3b1SBarry Smith 1350acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd() 135153acd3b1SBarry Smith 1352db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1353db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1354db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1355*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1356db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1357db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 135888aa4217SBarry Smith M*/ 135988aa4217SBarry Smith 13604416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 136153acd3b1SBarry Smith { 13624416b707SBarry Smith PetscOptionItem amsopt; 136353acd3b1SBarry Smith 136453acd3b1SBarry Smith PetscFunctionBegin; 1365e55864a3SBarry Smith if (!PetscOptionsObject->count) { 13669566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt)); 13679566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool),&amsopt->data)); 1368a297a907SKarl Rupp 1369ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 13701ae3d29cSBarry Smith } 137117326d04SJed Brown *flg = PETSC_FALSE; 13729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL)); 1373e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 13749566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man))); 137553acd3b1SBarry Smith } 137653acd3b1SBarry Smith PetscFunctionReturn(0); 137753acd3b1SBarry Smith } 137853acd3b1SBarry Smith 137988aa4217SBarry Smith /*MC 1380acfcf0e5SJed Brown PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 1381d5649816SBarry Smith which at most a single value can be true. 138253acd3b1SBarry Smith 13833f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 138453acd3b1SBarry Smith 138588aa4217SBarry Smith Synopsis: 138688aa4217SBarry Smith #include "petscsys.h" 13873a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool *flg) 138888aa4217SBarry Smith 138953acd3b1SBarry Smith Input Parameters: 139053acd3b1SBarry Smith + opt - option name 139153acd3b1SBarry Smith . text - short string that describes the option 139253acd3b1SBarry Smith - man - manual page with additional information on option 139353acd3b1SBarry Smith 139453acd3b1SBarry Smith Output Parameter: 139553acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 139653acd3b1SBarry Smith 139753acd3b1SBarry Smith Level: intermediate 139853acd3b1SBarry Smith 139995452b02SPatrick Sanan Notes: 140095452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 140153acd3b1SBarry Smith 1402acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() 140353acd3b1SBarry Smith 1404db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1405db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1406db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1407*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1408db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1409db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 141088aa4217SBarry Smith M*/ 141188aa4217SBarry Smith 14124416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 141353acd3b1SBarry Smith { 14144416b707SBarry Smith PetscOptionItem amsopt; 141553acd3b1SBarry Smith 141653acd3b1SBarry Smith PetscFunctionBegin; 1417e55864a3SBarry Smith if (!PetscOptionsObject->count) { 14189566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt)); 14199566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool),&amsopt->data)); 1420a297a907SKarl Rupp 1421ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 14221ae3d29cSBarry Smith } 142317326d04SJed Brown *flg = PETSC_FALSE; 14249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL)); 1425e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 14269566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man))); 142753acd3b1SBarry Smith } 142853acd3b1SBarry Smith PetscFunctionReturn(0); 142953acd3b1SBarry Smith } 143053acd3b1SBarry Smith 143188aa4217SBarry Smith /*MC 1432acfcf0e5SJed Brown PetscOptionsBool - Determines if a particular option is in the database with a true or false 143353acd3b1SBarry Smith 14343f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 143553acd3b1SBarry Smith 143688aa4217SBarry Smith Synopsis: 143788aa4217SBarry Smith #include "petscsys.h" 14383a89f35bSSatish Balay PetscErrorCode PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 143988aa4217SBarry Smith 144053acd3b1SBarry Smith Input Parameters: 144153acd3b1SBarry Smith + opt - option name 144253acd3b1SBarry Smith . text - short string that describes the option 1443868c398cSBarry Smith . man - manual page with additional information on option 144494ae4db5SBarry Smith - currentvalue - the current value 144553acd3b1SBarry Smith 1446d8d19677SJose E. Roman Output Parameters: 1447a2b725a8SWilliam Gropp + flg - PETSC_TRUE or PETSC_FALSE 1448a2b725a8SWilliam Gropp - set - PETSC_TRUE if found, else PETSC_FALSE 144953acd3b1SBarry Smith 14502efd9cb1SBarry Smith Notes: 14512efd9cb1SBarry Smith TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 14522efd9cb1SBarry Smith FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 14532efd9cb1SBarry Smith 14542efd9cb1SBarry 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 14552efd9cb1SBarry Smith is equivalent to -requested_bool true 14562efd9cb1SBarry Smith 14572efd9cb1SBarry Smith If the user does not supply the option at all flg is NOT changed. Thus 14582efd9cb1SBarry Smith you should ALWAYS initialize the flg if you access it without first checking if the set flag is true. 14592efd9cb1SBarry Smith 146095452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 146153acd3b1SBarry Smith 14625a856986SBarry Smith Level: beginner 14635a856986SBarry Smith 1464db781477SPatrick Sanan .seealso: `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1465db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsGetBool()`, 1466db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 1467db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1468*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1469db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1470db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 147188aa4217SBarry Smith M*/ 147288aa4217SBarry Smith 14734416b707SBarry Smith PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 147453acd3b1SBarry Smith { 1475ace3abfcSBarry Smith PetscBool iset; 14764416b707SBarry Smith PetscOptionItem amsopt; 147753acd3b1SBarry Smith 147853acd3b1SBarry Smith PetscFunctionBegin; 1479e55864a3SBarry Smith if (!PetscOptionsObject->count) { 14809566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt)); 14819566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(PetscBool),&amsopt->data)); 1482a297a907SKarl Rupp 148394ae4db5SBarry Smith *(PetscBool*)amsopt->data = currentvalue; 1484af6d86caSBarry Smith } 14859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset)); 148653acd3b1SBarry Smith if (set) *set = iset; 14871a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 148844ef3d73SBarry Smith const char *v = PetscBools[currentvalue], *vn = PetscBools[iset && flg ? *flg : currentvalue]; 14899566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: <%s : %s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,vn,text,ManSection(man))); 149053acd3b1SBarry Smith } 149153acd3b1SBarry Smith PetscFunctionReturn(0); 149253acd3b1SBarry Smith } 149353acd3b1SBarry Smith 149488aa4217SBarry Smith /*MC 149553acd3b1SBarry Smith PetscOptionsRealArray - Gets an array of double values for a particular 149653acd3b1SBarry Smith option in the database. The values must be separated with commas with 149753acd3b1SBarry Smith no intervening spaces. 149853acd3b1SBarry Smith 14993f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 150053acd3b1SBarry Smith 150188aa4217SBarry Smith Synopsis: 150288aa4217SBarry Smith #include "petscsys.h" 15033a89f35bSSatish Balay PetscErrorCode PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 150488aa4217SBarry Smith 150553acd3b1SBarry Smith Input Parameters: 150653acd3b1SBarry Smith + opt - the option one is seeking 150753acd3b1SBarry Smith . text - short string describing option 150853acd3b1SBarry Smith . man - manual page for option 1509c89788bdSBarry Smith - n - maximum number of values that value has room for 151053acd3b1SBarry Smith 1511d8d19677SJose E. Roman Output Parameters: 151253acd3b1SBarry Smith + value - location to copy values 1513c89788bdSBarry Smith . n - actual number of values found 151453acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 151553acd3b1SBarry Smith 151653acd3b1SBarry Smith Level: beginner 151753acd3b1SBarry Smith 151853acd3b1SBarry Smith Notes: 151953acd3b1SBarry Smith The user should pass in an array of doubles 152053acd3b1SBarry Smith 152153acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 152253acd3b1SBarry Smith 1523db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1524db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1525db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1526*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1527db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1528db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 152988aa4217SBarry Smith M*/ 153088aa4217SBarry Smith 15314416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 153253acd3b1SBarry Smith { 153353acd3b1SBarry Smith PetscInt i; 15344416b707SBarry Smith PetscOptionItem amsopt; 153553acd3b1SBarry Smith 153653acd3b1SBarry Smith PetscFunctionBegin; 1537e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1538e26ddf31SBarry Smith PetscReal *vals; 1539e26ddf31SBarry Smith 15409566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt)); 15419566063dSJacob Faibussowitsch PetscCall(PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data)); 1542e26ddf31SBarry Smith vals = (PetscReal*)amsopt->data; 1543e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1544e26ddf31SBarry Smith amsopt->arraylength = *n; 1545e26ddf31SBarry Smith } 15469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set)); 1547e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 15489566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0])); 154953acd3b1SBarry Smith for (i=1; i<*n; i++) { 15509566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i])); 155153acd3b1SBarry Smith } 15529566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man))); 155353acd3b1SBarry Smith } 155453acd3b1SBarry Smith PetscFunctionReturn(0); 155553acd3b1SBarry Smith } 155653acd3b1SBarry Smith 155788aa4217SBarry Smith /*MC 1558050cccc3SHong Zhang PetscOptionsScalarArray - Gets an array of Scalar values for a particular 1559050cccc3SHong Zhang option in the database. The values must be separated with commas with 1560050cccc3SHong Zhang no intervening spaces. 1561050cccc3SHong Zhang 1562050cccc3SHong Zhang Logically Collective on the communicator passed in PetscOptionsBegin() 1563050cccc3SHong Zhang 156488aa4217SBarry Smith Synopsis: 156588aa4217SBarry Smith #include "petscsys.h" 15663a89f35bSSatish Balay PetscErrorCode PetscOptionsScalarArray(const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool *set) 156788aa4217SBarry Smith 1568050cccc3SHong Zhang Input Parameters: 1569050cccc3SHong Zhang + opt - the option one is seeking 1570050cccc3SHong Zhang . text - short string describing option 1571050cccc3SHong Zhang . man - manual page for option 1572c89788bdSBarry Smith - n - maximum number of values allowed in the value array 1573050cccc3SHong Zhang 1574d8d19677SJose E. Roman Output Parameters: 1575050cccc3SHong Zhang + value - location to copy values 1576c89788bdSBarry Smith . n - actual number of values found 1577050cccc3SHong Zhang - set - PETSC_TRUE if found, else PETSC_FALSE 1578050cccc3SHong Zhang 1579050cccc3SHong Zhang Level: beginner 1580050cccc3SHong Zhang 1581050cccc3SHong Zhang Notes: 1582050cccc3SHong Zhang The user should pass in an array of doubles 1583050cccc3SHong Zhang 1584050cccc3SHong Zhang Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1585050cccc3SHong Zhang 1586db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1587db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1588db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1589*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1590db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1591db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 159288aa4217SBarry Smith M*/ 159388aa4217SBarry Smith 15944416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool *set) 1595050cccc3SHong Zhang { 1596050cccc3SHong Zhang PetscInt i; 15974416b707SBarry Smith PetscOptionItem amsopt; 1598050cccc3SHong Zhang 1599050cccc3SHong Zhang PetscFunctionBegin; 1600050cccc3SHong Zhang if (!PetscOptionsObject->count) { 1601050cccc3SHong Zhang PetscScalar *vals; 1602050cccc3SHong Zhang 16039566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt)); 16049566063dSJacob Faibussowitsch PetscCall(PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data)); 1605050cccc3SHong Zhang vals = (PetscScalar*)amsopt->data; 1606050cccc3SHong Zhang for (i=0; i<*n; i++) vals[i] = value[i]; 1607050cccc3SHong Zhang amsopt->arraylength = *n; 1608050cccc3SHong Zhang } 16099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set)); 1610050cccc3SHong Zhang if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 16119566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g+%gi",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)PetscRealPart(value[0]),(double)PetscImaginaryPart(value[0]))); 1612050cccc3SHong Zhang for (i=1; i<*n; i++) { 16139566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]))); 1614050cccc3SHong Zhang } 16159566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man))); 1616050cccc3SHong Zhang } 1617050cccc3SHong Zhang PetscFunctionReturn(0); 1618050cccc3SHong Zhang } 161953acd3b1SBarry Smith 162088aa4217SBarry Smith /*MC 162153acd3b1SBarry Smith PetscOptionsIntArray - Gets an array of integers for a particular 1622b32a342fSShri Abhyankar option in the database. 162353acd3b1SBarry Smith 16243f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 162553acd3b1SBarry Smith 162688aa4217SBarry Smith Synopsis: 162788aa4217SBarry Smith #include "petscsys.h" 16283a89f35bSSatish Balay PetscErrorCode PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 162988aa4217SBarry Smith 163053acd3b1SBarry Smith Input Parameters: 163153acd3b1SBarry Smith + opt - the option one is seeking 163253acd3b1SBarry Smith . text - short string describing option 163353acd3b1SBarry Smith . man - manual page for option 1634f8a50e2bSBarry Smith - n - maximum number of values 163553acd3b1SBarry Smith 1636d8d19677SJose E. Roman Output Parameters: 163753acd3b1SBarry Smith + value - location to copy values 1638f8a50e2bSBarry Smith . n - actual number of values found 163953acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 164053acd3b1SBarry Smith 164153acd3b1SBarry Smith Level: beginner 164253acd3b1SBarry Smith 164353acd3b1SBarry Smith Notes: 1644b32a342fSShri Abhyankar The array can be passed as 1645bebe2cf6SSatish Balay a comma separated list: 0,1,2,3,4,5,6,7 16460fd488f5SShri Abhyankar a range (start-end+1): 0-8 16470fd488f5SShri Abhyankar a range with given increment (start-end+1:inc): 0-7:2 1648bebe2cf6SSatish Balay a combination of values and ranges separated by commas: 0,1-8,8-15:2 1649b32a342fSShri Abhyankar 1650b32a342fSShri Abhyankar There must be no intervening spaces between the values. 165153acd3b1SBarry Smith 165253acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 165353acd3b1SBarry Smith 1654db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1655db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1656db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1657*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1658db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1659db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()`, `PetscOptionsRealArray()` 166088aa4217SBarry Smith M*/ 166188aa4217SBarry Smith 16624416b707SBarry Smith PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 166353acd3b1SBarry Smith { 166453acd3b1SBarry Smith PetscInt i; 16654416b707SBarry Smith PetscOptionItem amsopt; 166653acd3b1SBarry Smith 166753acd3b1SBarry Smith PetscFunctionBegin; 1668e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1669e26ddf31SBarry Smith PetscInt *vals; 1670e26ddf31SBarry Smith 16719566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt)); 16729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n,(PetscInt**)&amsopt->data)); 1673e26ddf31SBarry Smith vals = (PetscInt*)amsopt->data; 1674e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1675e26ddf31SBarry Smith amsopt->arraylength = *n; 1676e26ddf31SBarry Smith } 16779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set)); 1678e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 16799566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%" PetscInt_FMT,PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0])); 168053acd3b1SBarry Smith for (i=1; i<*n; i++) { 16819566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,",%" PetscInt_FMT,value[i])); 168253acd3b1SBarry Smith } 16839566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man))); 168453acd3b1SBarry Smith } 168553acd3b1SBarry Smith PetscFunctionReturn(0); 168653acd3b1SBarry Smith } 168753acd3b1SBarry Smith 168888aa4217SBarry Smith /*MC 168953acd3b1SBarry Smith PetscOptionsStringArray - Gets an array of string values for a particular 169053acd3b1SBarry Smith option in the database. The values must be separated with commas with 169153acd3b1SBarry Smith no intervening spaces. 169253acd3b1SBarry Smith 16933f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 169453acd3b1SBarry Smith 169588aa4217SBarry Smith Synopsis: 169688aa4217SBarry Smith #include "petscsys.h" 16973a89f35bSSatish Balay PetscErrorCode PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 169888aa4217SBarry Smith 169953acd3b1SBarry Smith Input Parameters: 170053acd3b1SBarry Smith + opt - the option one is seeking 170153acd3b1SBarry Smith . text - short string describing option 170253acd3b1SBarry Smith . man - manual page for option 170353acd3b1SBarry Smith - nmax - maximum number of strings 170453acd3b1SBarry Smith 1705d8d19677SJose E. Roman Output Parameters: 170653acd3b1SBarry Smith + value - location to copy strings 170753acd3b1SBarry Smith . nmax - actual number of strings found 170853acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 170953acd3b1SBarry Smith 171053acd3b1SBarry Smith Level: beginner 171153acd3b1SBarry Smith 171253acd3b1SBarry Smith Notes: 171353acd3b1SBarry Smith The user should pass in an array of pointers to char, to hold all the 171453acd3b1SBarry Smith strings returned by this function. 171553acd3b1SBarry Smith 171653acd3b1SBarry Smith The user is responsible for deallocating the strings that are 171753acd3b1SBarry Smith returned. The Fortran interface for this routine is not supported. 171853acd3b1SBarry Smith 171953acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 172053acd3b1SBarry Smith 1721db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1722db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1723db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1724*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1725db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1726db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 172788aa4217SBarry Smith M*/ 172888aa4217SBarry Smith 17294416b707SBarry Smith PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 173053acd3b1SBarry Smith { 17314416b707SBarry Smith PetscOptionItem amsopt; 173253acd3b1SBarry Smith 173353acd3b1SBarry Smith PetscFunctionBegin; 1734e55864a3SBarry Smith if (!PetscOptionsObject->count) { 17359566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt)); 17369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nmax,(char**)&amsopt->data)); 1737a297a907SKarl Rupp 17381ae3d29cSBarry Smith amsopt->arraylength = *nmax; 17391ae3d29cSBarry Smith } 17409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set)); 1741e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 17429566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man))); 174353acd3b1SBarry Smith } 174453acd3b1SBarry Smith PetscFunctionReturn(0); 174553acd3b1SBarry Smith } 174653acd3b1SBarry Smith 174788aa4217SBarry Smith /*MC 1748acfcf0e5SJed Brown PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 1749e2446a98SMatthew Knepley option in the database. The values must be separated with commas with 1750e2446a98SMatthew Knepley no intervening spaces. 1751e2446a98SMatthew Knepley 17523f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 1753e2446a98SMatthew Knepley 175488aa4217SBarry Smith Synopsis: 175588aa4217SBarry Smith #include "petscsys.h" 17563a89f35bSSatish Balay PetscErrorCode PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 175788aa4217SBarry Smith 1758e2446a98SMatthew Knepley Input Parameters: 1759e2446a98SMatthew Knepley + opt - the option one is seeking 1760e2446a98SMatthew Knepley . text - short string describing option 1761e2446a98SMatthew Knepley . man - manual page for option 1762c89788bdSBarry Smith - n - maximum number of values allowed in the value array 1763e2446a98SMatthew Knepley 1764d8d19677SJose E. Roman Output Parameters: 1765e2446a98SMatthew Knepley + value - location to copy values 1766c89788bdSBarry Smith . n - actual number of values found 1767e2446a98SMatthew Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 1768e2446a98SMatthew Knepley 1769e2446a98SMatthew Knepley Level: beginner 1770e2446a98SMatthew Knepley 1771e2446a98SMatthew Knepley Notes: 1772e2446a98SMatthew Knepley The user should pass in an array of doubles 1773e2446a98SMatthew Knepley 1774e2446a98SMatthew Knepley Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1775e2446a98SMatthew Knepley 1776db781477SPatrick Sanan .seealso: `PetscOptionsGetInt()`, `PetscOptionsGetReal()`, 1777db781477SPatrick Sanan `PetscOptionsHasName()`, `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`, 1778db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1779*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1780db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1781db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 178288aa4217SBarry Smith M*/ 178388aa4217SBarry Smith 17844416b707SBarry Smith PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 1785e2446a98SMatthew Knepley { 1786e2446a98SMatthew Knepley PetscInt i; 17874416b707SBarry Smith PetscOptionItem amsopt; 1788e2446a98SMatthew Knepley 1789e2446a98SMatthew Knepley PetscFunctionBegin; 1790e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1791ace3abfcSBarry Smith PetscBool *vals; 17921ae3d29cSBarry Smith 17939566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt)); 17949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*n,(PetscBool**)&amsopt->data)); 1795ace3abfcSBarry Smith vals = (PetscBool*)amsopt->data; 17961ae3d29cSBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 17971ae3d29cSBarry Smith amsopt->arraylength = *n; 17981ae3d29cSBarry Smith } 17999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set)); 1800e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 18019566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0])); 1802e2446a98SMatthew Knepley for (i=1; i<*n; i++) { 18039566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i])); 1804e2446a98SMatthew Knepley } 18059566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man))); 1806e2446a98SMatthew Knepley } 1807e2446a98SMatthew Knepley PetscFunctionReturn(0); 1808e2446a98SMatthew Knepley } 1809e2446a98SMatthew Knepley 181088aa4217SBarry Smith /*MC 1811d1da0b69SBarry Smith PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user 18128cc676e6SMatthew G Knepley 18138cc676e6SMatthew G Knepley Logically Collective on the communicator passed in PetscOptionsBegin() 18148cc676e6SMatthew G Knepley 181588aa4217SBarry Smith Synopsis: 181688aa4217SBarry Smith #include "petscsys.h" 18173a89f35bSSatish Balay PetscErrorCode PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 181888aa4217SBarry Smith 18198cc676e6SMatthew G Knepley Input Parameters: 18208cc676e6SMatthew G Knepley + opt - option name 18218cc676e6SMatthew G Knepley . text - short string that describes the option 18228cc676e6SMatthew G Knepley - man - manual page with additional information on option 18238cc676e6SMatthew G Knepley 1824d8d19677SJose E. Roman Output Parameters: 18258cc676e6SMatthew G Knepley + viewer - the viewer 18267962402dSFande Kong . format - the PetscViewerFormat requested by the user, pass NULL if not needed 18278cc676e6SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 18288cc676e6SMatthew G Knepley 18298cc676e6SMatthew G Knepley Level: beginner 18308cc676e6SMatthew G Knepley 183195452b02SPatrick Sanan Notes: 183295452b02SPatrick Sanan Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 18338cc676e6SMatthew G Knepley 18345a7113b9SPatrick Sanan See PetscOptionsGetViewer() for the format of the supplied viewer and its options 18358cc676e6SMatthew G Knepley 1836db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, `PetscOptionsGetInt()`, 1837db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1838db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1839db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1840*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1841db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1842db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 184388aa4217SBarry Smith M*/ 184488aa4217SBarry Smith 18454416b707SBarry Smith PetscErrorCode PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 18468cc676e6SMatthew G Knepley { 18474416b707SBarry Smith PetscOptionItem amsopt; 18488cc676e6SMatthew G Knepley 18498cc676e6SMatthew G Knepley PetscFunctionBegin; 18501a1499c8SBarry Smith if (!PetscOptionsObject->count) { 18519566063dSJacob Faibussowitsch PetscCall(PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt)); 185264facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 18539566063dSJacob Faibussowitsch PetscCall(PetscStrdup("",(char**)&amsopt->data)); 18548cc676e6SMatthew G Knepley } 18559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->options,PetscOptionsObject->prefix,opt,viewer,format,set)); 1856e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 18579566063dSJacob Faibussowitsch PetscCall((*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man))); 18588cc676e6SMatthew G Knepley } 18598cc676e6SMatthew G Knepley PetscFunctionReturn(0); 18608cc676e6SMatthew G Knepley } 1861