17d0a6c19SBarry Smith 21a1499c8SBarry Smith 353acd3b1SBarry Smith /* 43fc1eb6aSBarry Smith Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to 53fc1eb6aSBarry Smith GUI code to display the options and get values from the users. 653acd3b1SBarry Smith 753acd3b1SBarry Smith */ 853acd3b1SBarry Smith 9afcb2eb5SJed Brown #include <petsc-private/petscimpl.h> /*I "petscsys.h" I*/ 10665c2dedSJed Brown #include <petscviewer.h> 1153acd3b1SBarry Smith 122aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None") 132aa6d131SJed Brown 1453acd3b1SBarry Smith /* 1553acd3b1SBarry Smith Keep a linked list of options that have been posted and we are waiting for 163fc1eb6aSBarry Smith user selection. See the manual page for PetscOptionsBegin() 1753acd3b1SBarry Smith 1853acd3b1SBarry Smith Eventually we'll attach this beast to a MPI_Comm 1953acd3b1SBarry Smith */ 20e55864a3SBarry Smith 2153acd3b1SBarry Smith 2253acd3b1SBarry Smith #undef __FUNCT__ 2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private" 2453acd3b1SBarry Smith /* 2553acd3b1SBarry Smith Handles setting up the data structure in a call to PetscOptionsBegin() 2653acd3b1SBarry Smith */ 278c34d3f5SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptions *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[]) 2853acd3b1SBarry Smith { 2953acd3b1SBarry Smith PetscErrorCode ierr; 3053acd3b1SBarry Smith 3153acd3b1SBarry Smith PetscFunctionBegin; 32e55864a3SBarry Smith PetscOptionsObject->next = 0; 33e55864a3SBarry Smith PetscOptionsObject->comm = comm; 34e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_FALSE; 35a297a907SKarl Rupp 36e55864a3SBarry Smith ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr); 37e55864a3SBarry Smith ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr); 3853acd3b1SBarry Smith 39e55864a3SBarry Smith ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr); 40e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) { 41e55864a3SBarry Smith if (!PetscOptionsObject->alreadyprinted) { 4253acd3b1SBarry Smith ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr); 4353acd3b1SBarry Smith } 4461b37b28SSatish Balay } 4553acd3b1SBarry Smith PetscFunctionReturn(0); 4653acd3b1SBarry Smith } 4753acd3b1SBarry Smith 483194b578SJed Brown #undef __FUNCT__ 493194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private" 503194b578SJed Brown /* 513194b578SJed Brown Handles setting up the data structure in a call to PetscObjectOptionsBegin() 523194b578SJed Brown */ 538c34d3f5SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptions *PetscOptionsObject,PetscObject obj) 543194b578SJed Brown { 553194b578SJed Brown PetscErrorCode ierr; 563194b578SJed Brown char title[256]; 573194b578SJed Brown PetscBool flg; 583194b578SJed Brown 593194b578SJed Brown PetscFunctionBegin; 603194b578SJed Brown PetscValidHeader(obj,1); 61e55864a3SBarry Smith PetscOptionsObject->object = obj; 62e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = obj->optionsprinted; 63a297a907SKarl Rupp 643194b578SJed Brown ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr); 653194b578SJed Brown if (flg) { 668caf3d72SBarry Smith ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr); 673194b578SJed Brown } else { 688caf3d72SBarry Smith ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr); 693194b578SJed Brown } 70e55864a3SBarry Smith ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr); 713194b578SJed Brown PetscFunctionReturn(0); 723194b578SJed Brown } 733194b578SJed Brown 7453acd3b1SBarry Smith /* 7553acd3b1SBarry Smith Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd() 7653acd3b1SBarry Smith */ 7753acd3b1SBarry Smith #undef __FUNCT__ 7853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private" 798c34d3f5SBarry Smith static int PetscOptionsCreate_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOption *amsopt) 8053acd3b1SBarry Smith { 8153acd3b1SBarry Smith int ierr; 828c34d3f5SBarry Smith PetscOption next; 833be6e4c3SJed Brown PetscBool valid; 8453acd3b1SBarry Smith 8553acd3b1SBarry Smith PetscFunctionBegin; 863be6e4c3SJed Brown ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr); 873be6e4c3SJed Brown if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt); 883be6e4c3SJed Brown 89b00a9115SJed Brown ierr = PetscNew(amsopt);CHKERRQ(ierr); 9053acd3b1SBarry Smith (*amsopt)->next = 0; 9153acd3b1SBarry Smith (*amsopt)->set = PETSC_FALSE; 926356e834SBarry Smith (*amsopt)->type = t; 9353acd3b1SBarry Smith (*amsopt)->data = 0; 9461b37b28SSatish Balay 9553acd3b1SBarry Smith ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr); 9653acd3b1SBarry Smith ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr); 976356e834SBarry Smith ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr); 9853acd3b1SBarry Smith 99e55864a3SBarry Smith if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt; 100a297a907SKarl Rupp else { 101e55864a3SBarry Smith next = PetscOptionsObject->next; 10253acd3b1SBarry Smith while (next->next) next = next->next; 10353acd3b1SBarry Smith next->next = *amsopt; 10453acd3b1SBarry Smith } 10553acd3b1SBarry Smith PetscFunctionReturn(0); 10653acd3b1SBarry Smith } 10753acd3b1SBarry Smith 10853acd3b1SBarry Smith #undef __FUNCT__ 109aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString" 110aee2cecaSBarry Smith /* 1113fc1eb6aSBarry Smith PetscScanString - Gets user input via stdin from process and broadcasts to all processes 1123fc1eb6aSBarry Smith 1133fc1eb6aSBarry Smith Collective on MPI_Comm 1143fc1eb6aSBarry Smith 1153fc1eb6aSBarry Smith Input Parameters: 1163fc1eb6aSBarry Smith + commm - communicator for the broadcast, must be PETSC_COMM_WORLD 1173fc1eb6aSBarry Smith . n - length of the string, must be the same on all processes 1183fc1eb6aSBarry Smith - str - location to store input 119aee2cecaSBarry Smith 120aee2cecaSBarry Smith Bugs: 121aee2cecaSBarry Smith . Assumes process 0 of the given communicator has access to stdin 122aee2cecaSBarry Smith 123aee2cecaSBarry Smith */ 1243fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[]) 125aee2cecaSBarry Smith { 126330cf3c9SBarry Smith size_t i; 127aee2cecaSBarry Smith char c; 1283fc1eb6aSBarry Smith PetscMPIInt rank,nm; 129aee2cecaSBarry Smith PetscErrorCode ierr; 130aee2cecaSBarry Smith 131aee2cecaSBarry Smith PetscFunctionBegin; 132aee2cecaSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 133aee2cecaSBarry Smith if (!rank) { 134aee2cecaSBarry Smith c = (char) getchar(); 135aee2cecaSBarry Smith i = 0; 136aee2cecaSBarry Smith while (c != '\n' && i < n-1) { 137aee2cecaSBarry Smith str[i++] = c; 138aee2cecaSBarry Smith c = (char)getchar(); 139aee2cecaSBarry Smith } 140aee2cecaSBarry Smith str[i] = 0; 141aee2cecaSBarry Smith } 1424dc2109aSBarry Smith ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr); 1433fc1eb6aSBarry Smith ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr); 144aee2cecaSBarry Smith PetscFunctionReturn(0); 145aee2cecaSBarry Smith } 146aee2cecaSBarry Smith 147ead66b60SBarry Smith #undef __FUNCT__ 148ead66b60SBarry Smith #define __FUNCT__ "PetscStrdup" 1495b02f95dSBarry Smith /* 1505b02f95dSBarry Smith This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy() 1515b02f95dSBarry Smith */ 1525b02f95dSBarry Smith static PetscErrorCode PetscStrdup(const char s[],char *t[]) 1535b02f95dSBarry Smith { 1545b02f95dSBarry Smith PetscErrorCode ierr; 1555b02f95dSBarry Smith size_t len; 1565b02f95dSBarry Smith char *tmp = 0; 1575b02f95dSBarry Smith 1585b02f95dSBarry Smith PetscFunctionBegin; 1595b02f95dSBarry Smith if (s) { 1605b02f95dSBarry Smith ierr = PetscStrlen(s,&len);CHKERRQ(ierr); 161ee48884fSBarry Smith tmp = (char*) malloc((len+1)*sizeof(char*)); 1625b02f95dSBarry Smith if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string"); 1635b02f95dSBarry Smith ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr); 1645b02f95dSBarry Smith } 1655b02f95dSBarry Smith *t = tmp; 1665b02f95dSBarry Smith PetscFunctionReturn(0); 1675b02f95dSBarry Smith } 1685b02f95dSBarry Smith 1695b02f95dSBarry Smith 170aee2cecaSBarry Smith #undef __FUNCT__ 171aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput" 172aee2cecaSBarry Smith /* 1733cc1e11dSBarry Smith PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime 174aee2cecaSBarry Smith 175aee2cecaSBarry Smith Notes: this isn't really practical, it is just to demonstrate the principle 176aee2cecaSBarry Smith 1777781c08eSBarry Smith A carriage return indicates no change from the default; but this like -ksp_monitor <stdout> the default is actually not stdout the default 1787781c08eSBarry Smith is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug? 1797781c08eSBarry Smith 180aee2cecaSBarry Smith Bugs: 1817781c08eSBarry Smith + All processes must traverse through the exact same set of option queries due to the call to PetscScanString() 1823cc1e11dSBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 183aee2cecaSBarry Smith - Only works for PetscInt == int, PetscReal == double etc 184aee2cecaSBarry Smith 1853cc1e11dSBarry Smith Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different 1863cc1e11dSBarry Smith address space and communicating with the PETSc program 1873cc1e11dSBarry Smith 188aee2cecaSBarry Smith */ 1898c34d3f5SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptions *PetscOptionsObject) 1906356e834SBarry Smith { 1916356e834SBarry Smith PetscErrorCode ierr; 1928c34d3f5SBarry Smith PetscOption next = PetscOptionsObject->next; 1936356e834SBarry Smith char str[512]; 1947781c08eSBarry Smith PetscBool bid; 195a4404d99SBarry Smith PetscReal ir,*valr; 196330cf3c9SBarry Smith PetscInt *vald; 197330cf3c9SBarry Smith size_t i; 1986356e834SBarry Smith 199e55864a3SBarry Smith ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr); 2006356e834SBarry Smith while (next) { 2016356e834SBarry Smith switch (next->type) { 2026356e834SBarry Smith case OPTION_HEAD: 2036356e834SBarry Smith break; 204e26ddf31SBarry Smith case OPTION_INT_ARRAY: 205e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 206e26ddf31SBarry Smith vald = (PetscInt*) next->data; 207e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 208e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr); 209e26ddf31SBarry Smith if (i < next->arraylength-1) { 210e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 211e26ddf31SBarry Smith } 212e26ddf31SBarry Smith } 213e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 214e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 215e26ddf31SBarry Smith if (str[0]) { 216e26ddf31SBarry Smith PetscToken token; 217e26ddf31SBarry Smith PetscInt n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end; 218e26ddf31SBarry Smith size_t len; 219e26ddf31SBarry Smith char *value; 220ace3abfcSBarry Smith PetscBool foundrange; 221e26ddf31SBarry Smith 222e26ddf31SBarry Smith next->set = PETSC_TRUE; 223e26ddf31SBarry Smith value = str; 224e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 225e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 226e26ddf31SBarry Smith while (n < nmax) { 227e26ddf31SBarry Smith if (!value) break; 228e26ddf31SBarry Smith 229e26ddf31SBarry Smith /* look for form d-D where d and D are integers */ 230e26ddf31SBarry Smith foundrange = PETSC_FALSE; 231e26ddf31SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 232e26ddf31SBarry Smith if (value[0] == '-') i=2; 233e26ddf31SBarry Smith else i=1; 234330cf3c9SBarry Smith for (;i<len; i++) { 235e26ddf31SBarry Smith if (value[i] == '-') { 236e32f2f54SBarry Smith if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 237e26ddf31SBarry Smith value[i] = 0; 238cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 239cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 240e32f2f54SBarry Smith if (end <= start) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, %s-%s cannot have decreasing list",n,value,value+i+1); 241e32f2f54SBarry Smith if (n + end - start - 1 >= nmax) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, not enough space in left in array (%D) to contain entire range from %D to %D",n,nmax-n,start,end); 242e26ddf31SBarry Smith for (; start<end; start++) { 243e26ddf31SBarry Smith *dvalue = start; dvalue++;n++; 244e26ddf31SBarry Smith } 245e26ddf31SBarry Smith foundrange = PETSC_TRUE; 246e26ddf31SBarry Smith break; 247e26ddf31SBarry Smith } 248e26ddf31SBarry Smith } 249e26ddf31SBarry Smith if (!foundrange) { 250cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr); 251e26ddf31SBarry Smith dvalue++; 252e26ddf31SBarry Smith n++; 253e26ddf31SBarry Smith } 254e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 255e26ddf31SBarry Smith } 2568c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 257e26ddf31SBarry Smith } 258e26ddf31SBarry Smith break; 259e26ddf31SBarry Smith case OPTION_REAL_ARRAY: 260e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 261e26ddf31SBarry Smith valr = (PetscReal*) next->data; 262e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 263e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr); 264e26ddf31SBarry Smith if (i < next->arraylength-1) { 265e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 266e26ddf31SBarry Smith } 267e26ddf31SBarry Smith } 268e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 269e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 270e26ddf31SBarry Smith if (str[0]) { 271e26ddf31SBarry Smith PetscToken token; 272e26ddf31SBarry Smith PetscInt n = 0,nmax = next->arraylength; 273e26ddf31SBarry Smith PetscReal *dvalue = (PetscReal*)next->data; 274e26ddf31SBarry Smith char *value; 275e26ddf31SBarry Smith 276e26ddf31SBarry Smith next->set = PETSC_TRUE; 277e26ddf31SBarry Smith value = str; 278e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 279e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 280e26ddf31SBarry Smith while (n < nmax) { 281e26ddf31SBarry Smith if (!value) break; 282cfbddea1SSatish Balay ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 283e26ddf31SBarry Smith dvalue++; 284e26ddf31SBarry Smith n++; 285e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 286e26ddf31SBarry Smith } 2878c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 288e26ddf31SBarry Smith } 289e26ddf31SBarry Smith break; 2906356e834SBarry Smith case OPTION_INT: 291e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%d>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,*(int*)next->data,next->text,next->man);CHKERRQ(ierr); 2923fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 2933fc1eb6aSBarry Smith if (str[0]) { 294d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 295d25d7f95SJed Brown long long lid; 296d25d7f95SJed Brown sscanf(str,"%lld",&lid); 2971a1499c8SBarry Smith if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %lld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid); 298c272547aSJed Brown #else 299d25d7f95SJed Brown long lid; 300d25d7f95SJed Brown sscanf(str,"%ld",&lid); 3011a1499c8SBarry Smith if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %ld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid); 302c272547aSJed Brown #endif 303a297a907SKarl Rupp 304d25d7f95SJed Brown next->set = PETSC_TRUE; 305d25d7f95SJed Brown *((PetscInt*)next->data) = (PetscInt)lid; 306aee2cecaSBarry Smith } 307aee2cecaSBarry Smith break; 308aee2cecaSBarry Smith case OPTION_REAL: 309e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%g>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,*(double*)next->data,next->text,next->man);CHKERRQ(ierr); 3103fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3113fc1eb6aSBarry Smith if (str[0]) { 312ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 313a4404d99SBarry Smith sscanf(str,"%e",&ir); 314ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE) 315aee2cecaSBarry Smith sscanf(str,"%le",&ir); 316ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 317d9822059SBarry Smith ir = strtoflt128(str,0); 318d9822059SBarry Smith #else 319513dbe71SLisandro Dalcin SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type"); 320a4404d99SBarry Smith #endif 321aee2cecaSBarry Smith next->set = PETSC_TRUE; 322aee2cecaSBarry Smith *((PetscReal*)next->data) = ir; 323aee2cecaSBarry Smith } 324aee2cecaSBarry Smith break; 3257781c08eSBarry Smith case OPTION_BOOL: 32683355fc5SBarry Smith ierr = 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);CHKERRQ(ierr); 3277781c08eSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3287781c08eSBarry Smith if (str[0]) { 3297781c08eSBarry Smith ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr); 3307781c08eSBarry Smith next->set = PETSC_TRUE; 3317781c08eSBarry Smith *((PetscBool*)next->data) = bid; 3327781c08eSBarry Smith } 3337781c08eSBarry Smith break; 334aee2cecaSBarry Smith case OPTION_STRING: 335e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,(char*)next->data,next->text,next->man);CHKERRQ(ierr); 3363fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3373fc1eb6aSBarry Smith if (str[0]) { 338aee2cecaSBarry Smith next->set = PETSC_TRUE; 33964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3405b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3416356e834SBarry Smith } 3426356e834SBarry Smith break; 343a264d7a6SBarry Smith case OPTION_FLIST: 344e55864a3SBarry Smith ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr); 3453cc1e11dSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3463cc1e11dSBarry Smith if (str[0]) { 347e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_TRUE; 3483cc1e11dSBarry Smith next->set = PETSC_TRUE; 34964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3505b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3513cc1e11dSBarry Smith } 3523cc1e11dSBarry Smith break; 353b432afdaSMatthew Knepley default: 354b432afdaSMatthew Knepley break; 3556356e834SBarry Smith } 3566356e834SBarry Smith next = next->next; 3576356e834SBarry Smith } 3586356e834SBarry Smith PetscFunctionReturn(0); 3596356e834SBarry Smith } 3606356e834SBarry Smith 361e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 362e04113cfSBarry Smith #include <petscviewersaws.h> 363d5649816SBarry Smith 364d5649816SBarry Smith static int count = 0; 365d5649816SBarry Smith 366b3506946SBarry Smith #undef __FUNCT__ 367e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy" 368e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void) 369d5649816SBarry Smith { 3702657e9d9SBarry Smith PetscFunctionBegin; 371d5649816SBarry Smith PetscFunctionReturn(0); 372d5649816SBarry Smith } 373d5649816SBarry Smith 3749c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n" 37523a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n" 37623a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n" 377d1fc0251SBarry Smith "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n" 37864bbc9efSBarry Smith "<script>\n" 37964bbc9efSBarry Smith "jQuery(document).ready(function() {\n" 38064bbc9efSBarry Smith "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n" 38164bbc9efSBarry Smith "})\n" 38264bbc9efSBarry Smith "</script>\n" 38364bbc9efSBarry Smith "</head>\n"; 384*1423471aSBarry Smith 385*1423471aSBarry Smith /* Determines the size and style of the scroll region where PETSc options selectable from users are displayed */ 386*1423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>"; 38764bbc9efSBarry Smith 388d5649816SBarry Smith #undef __FUNCT__ 3897781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput" 390b3506946SBarry Smith /* 3917781c08eSBarry Smith PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs 392b3506946SBarry Smith 393b3506946SBarry Smith Bugs: 394b3506946SBarry Smith + All processes must traverse through the exact same set of option queries do to the call to PetscScanString() 395b3506946SBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 396b3506946SBarry Smith - Only works for PetscInt == int, PetscReal == double etc 397b3506946SBarry Smith 398b3506946SBarry Smith 399b3506946SBarry Smith */ 400f7b25cf6SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptions *PetscOptionsObject) 401b3506946SBarry Smith { 402b3506946SBarry Smith PetscErrorCode ierr; 4038c34d3f5SBarry Smith PetscOption next = PetscOptionsObject->next; 404d5649816SBarry Smith static int mancount = 0; 405b3506946SBarry Smith char options[16]; 4067aab2a10SBarry Smith PetscBool changedmethod = PETSC_FALSE; 407a23eb982SSurtai Han PetscBool stopasking = PETSC_FALSE; 40888a9752cSBarry Smith char manname[16],textname[16]; 4092657e9d9SBarry Smith char dir[1024]; 410b3506946SBarry Smith 411b3506946SBarry Smith /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */ 412b3506946SBarry Smith sprintf(options,"Options_%d",count++); 413a297a907SKarl Rupp 414e55864a3SBarry Smith PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* AMS will change this, so cannot pass prefix directly */ 4151bc75a8dSBarry Smith 416d50831c4SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr); 41783355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING)); 4187781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr); 41983355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING)); 4202657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN)); 421a23eb982SSurtai Han PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN)); 422b3506946SBarry Smith 423b3506946SBarry Smith while (next) { 424d50831c4SBarry Smith sprintf(manname,"_man_%d",mancount); 4252657e9d9SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr); 4267781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING)); 427d50831c4SBarry Smith sprintf(textname,"_text_%d",mancount++); 4287781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr); 4297781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING)); 4309f32e415SBarry Smith 431b3506946SBarry Smith switch (next->type) { 432b3506946SBarry Smith case OPTION_HEAD: 433b3506946SBarry Smith break; 434b3506946SBarry Smith case OPTION_INT_ARRAY: 4357781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4362657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT)); 437b3506946SBarry Smith break; 438b3506946SBarry Smith case OPTION_REAL_ARRAY: 4397781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4402657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE)); 441b3506946SBarry Smith break; 442b3506946SBarry Smith case OPTION_INT: 4437781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4442657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT)); 445b3506946SBarry Smith break; 446b3506946SBarry Smith case OPTION_REAL: 4477781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4482657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE)); 449b3506946SBarry Smith break; 4507781c08eSBarry Smith case OPTION_BOOL: 4517781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4522657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN)); 4531ae3d29cSBarry Smith break; 4547781c08eSBarry Smith case OPTION_BOOL_ARRAY: 4557781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4562657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN)); 45771f08665SBarry Smith break; 458b3506946SBarry Smith case OPTION_STRING: 4597781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4607781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4611ae3d29cSBarry Smith break; 4621ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 4637781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4642657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING)); 465b3506946SBarry Smith break; 466a264d7a6SBarry Smith case OPTION_FLIST: 467a264d7a6SBarry Smith { 468a264d7a6SBarry Smith PetscInt ntext; 4697781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4707781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 471a264d7a6SBarry Smith ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr); 472a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 473a264d7a6SBarry Smith } 474a264d7a6SBarry Smith break; 4751ae3d29cSBarry Smith case OPTION_ELIST: 476a264d7a6SBarry Smith { 477a264d7a6SBarry Smith PetscInt ntext = next->nlist; 4787781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4797781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 480ead66b60SBarry Smith ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr); 4811ae3d29cSBarry Smith ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr); 482a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 483a264d7a6SBarry Smith } 484a264d7a6SBarry Smith break; 485b3506946SBarry Smith default: 486b3506946SBarry Smith break; 487b3506946SBarry Smith } 488b3506946SBarry Smith next = next->next; 489b3506946SBarry Smith } 490b3506946SBarry Smith 491b3506946SBarry Smith /* wait until accessor has unlocked the memory */ 49264bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader)); 49364bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom)); 4947aab2a10SBarry Smith ierr = PetscSAWsBlock();CHKERRQ(ierr); 49564bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Header,("index.html")); 49664bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2)); 497b3506946SBarry Smith 49888a9752cSBarry Smith /* determine if any values have been set in GUI */ 49983355fc5SBarry Smith next = PetscOptionsObject->next; 50088a9752cSBarry Smith while (next) { 50188a9752cSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 502f7b25cf6SBarry Smith PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set)); 50388a9752cSBarry Smith next = next->next; 50488a9752cSBarry Smith } 50588a9752cSBarry Smith 506b3506946SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 507f7b25cf6SBarry Smith if (changedmethod) PetscOptionsObject->count = -2; 508b3506946SBarry Smith 509a23eb982SSurtai Han if (stopasking) { 510a23eb982SSurtai Han PetscOptionsPublish = PETSC_FALSE; 51112655325SBarry Smith PetscOptionsObject->count = 0;//do not ask for same thing again 512a23eb982SSurtai Han } 513a23eb982SSurtai Han 5149a492a5cSBarry Smith PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options")); 515b3506946SBarry Smith PetscFunctionReturn(0); 516b3506946SBarry Smith } 517b3506946SBarry Smith #endif 518b3506946SBarry Smith 5196356e834SBarry Smith #undef __FUNCT__ 52053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private" 5218c34d3f5SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptions *PetscOptionsObject) 52253acd3b1SBarry Smith { 52353acd3b1SBarry Smith PetscErrorCode ierr; 5248c34d3f5SBarry Smith PetscOption last; 5256356e834SBarry Smith char option[256],value[1024],tmp[32]; 526330cf3c9SBarry Smith size_t j; 52753acd3b1SBarry Smith 52853acd3b1SBarry Smith PetscFunctionBegin; 52983355fc5SBarry Smith if (PetscOptionsObject->next) { 53083355fc5SBarry Smith if (!PetscOptionsObject->count) { 531a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS) 532f7b25cf6SBarry Smith ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr); 533b3506946SBarry Smith #else 534e55864a3SBarry Smith ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr); 535b3506946SBarry Smith #endif 536aee2cecaSBarry Smith } 537aee2cecaSBarry Smith } 5386356e834SBarry Smith 539e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr); 540e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr); 5416356e834SBarry Smith 542e26ddf31SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 543e55864a3SBarry Smith if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2; 5447a72a596SBarry Smith /* reset alreadyprinted flag */ 545e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = PETSC_FALSE; 546e55864a3SBarry Smith if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE; 547e55864a3SBarry Smith PetscOptionsObject->object = NULL; 54853acd3b1SBarry Smith 549e55864a3SBarry Smith while (PetscOptionsObject->next) { 550e55864a3SBarry Smith if (PetscOptionsObject->next->set) { 551e55864a3SBarry Smith if (PetscOptionsObject->prefix) { 55253acd3b1SBarry Smith ierr = PetscStrcpy(option,"-");CHKERRQ(ierr); 553e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr); 554e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr); 5556356e834SBarry Smith } else { 556e55864a3SBarry Smith ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr); 5576356e834SBarry Smith } 5586356e834SBarry Smith 559e55864a3SBarry Smith switch (PetscOptionsObject->next->type) { 5606356e834SBarry Smith case OPTION_HEAD: 5616356e834SBarry Smith break; 5626356e834SBarry Smith case OPTION_INT_ARRAY: 563e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]); 564e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 565e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]); 5666356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5676356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5686356e834SBarry Smith } 5696356e834SBarry Smith break; 5706356e834SBarry Smith case OPTION_INT: 571e55864a3SBarry Smith sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data); 5726356e834SBarry Smith break; 5736356e834SBarry Smith case OPTION_REAL: 574e55864a3SBarry Smith sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data); 5756356e834SBarry Smith break; 5766356e834SBarry Smith case OPTION_REAL_ARRAY: 577e55864a3SBarry Smith sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]); 578e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 579e55864a3SBarry Smith sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]); 5806356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5816356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5826356e834SBarry Smith } 5836356e834SBarry Smith break; 5847781c08eSBarry Smith case OPTION_BOOL: 585e55864a3SBarry Smith sprintf(value,"%d",*(int*)PetscOptionsObject->next->data); 5866356e834SBarry Smith break; 5877781c08eSBarry Smith case OPTION_BOOL_ARRAY: 588e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]); 589e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 590e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]); 5911ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5921ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5931ae3d29cSBarry Smith } 5941ae3d29cSBarry Smith break; 595a264d7a6SBarry Smith case OPTION_FLIST: 5961ae3d29cSBarry Smith case OPTION_ELIST: 597e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 5986356e834SBarry Smith break; 5991ae3d29cSBarry Smith case OPTION_STRING: 600e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 6011ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 602e55864a3SBarry Smith sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]); 603e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 604e55864a3SBarry Smith sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]); 6051ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 6061ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 6071ae3d29cSBarry Smith } 6086356e834SBarry Smith break; 6096356e834SBarry Smith } 6106356e834SBarry Smith ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr); 6116356e834SBarry Smith } 612e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr); 613e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr); 614e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr); 615e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr); 616c979a496SBarry Smith 61783355fc5SBarry Smith if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){ 61883355fc5SBarry Smith free(PetscOptionsObject->next->data); 619c979a496SBarry Smith } else { 62083355fc5SBarry Smith ierr = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr); 621c979a496SBarry Smith } 6227781c08eSBarry Smith 62383355fc5SBarry Smith last = PetscOptionsObject->next; 62483355fc5SBarry Smith PetscOptionsObject->next = PetscOptionsObject->next->next; 6256356e834SBarry Smith ierr = PetscFree(last);CHKERRQ(ierr); 6266356e834SBarry Smith } 627e55864a3SBarry Smith PetscOptionsObject->next = 0; 62853acd3b1SBarry Smith PetscFunctionReturn(0); 62953acd3b1SBarry Smith } 63053acd3b1SBarry Smith 63153acd3b1SBarry Smith #undef __FUNCT__ 632e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private" 63353acd3b1SBarry Smith /*@C 63453acd3b1SBarry Smith PetscOptionsEnum - Gets the enum value for a particular option in the database. 63553acd3b1SBarry Smith 6363f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 63753acd3b1SBarry Smith 63853acd3b1SBarry Smith Input Parameters: 63953acd3b1SBarry Smith + opt - option name 64053acd3b1SBarry Smith . text - short string that describes the option 64153acd3b1SBarry Smith . man - manual page with additional information on option 64253acd3b1SBarry Smith . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 6430fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6440fdccdaeSBarry Smith $ PetscOptionsEnum(..., obj->value,&object->value,...) or 6450fdccdaeSBarry Smith $ value = defaultvalue 6460fdccdaeSBarry Smith $ PetscOptionsEnum(..., value,&value,&flg); 6470fdccdaeSBarry Smith $ if (flg) { 64853acd3b1SBarry Smith 64953acd3b1SBarry Smith Output Parameter: 65053acd3b1SBarry Smith + value - the value to return 651b32e0204SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 65253acd3b1SBarry Smith 65353acd3b1SBarry Smith Level: beginner 65453acd3b1SBarry Smith 65553acd3b1SBarry Smith Concepts: options database 65653acd3b1SBarry Smith 65753acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 65853acd3b1SBarry Smith 65953acd3b1SBarry Smith list is usually something like PCASMTypes or some other predefined list of enum names 66053acd3b1SBarry Smith 66153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 662acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 663acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 66453acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 66553acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 666acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 667a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 66853acd3b1SBarry Smith @*/ 6698c34d3f5SBarry Smith PetscErrorCode PetscOptionsEnum_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool *set) 67053acd3b1SBarry Smith { 67153acd3b1SBarry Smith PetscErrorCode ierr; 67253acd3b1SBarry Smith PetscInt ntext = 0; 673aa5bb8c0SSatish Balay PetscInt tval; 674ace3abfcSBarry Smith PetscBool tflg; 67553acd3b1SBarry Smith 67653acd3b1SBarry Smith PetscFunctionBegin; 67753acd3b1SBarry Smith while (list[ntext++]) { 678e32f2f54SBarry Smith if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 67953acd3b1SBarry Smith } 680e32f2f54SBarry Smith if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 68153acd3b1SBarry Smith ntext -= 3; 682e55864a3SBarry Smith ierr = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr); 683aa5bb8c0SSatish Balay /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 684aa5bb8c0SSatish Balay if (tflg) *value = (PetscEnum)tval; 685aa5bb8c0SSatish Balay if (set) *set = tflg; 68653acd3b1SBarry Smith PetscFunctionReturn(0); 68753acd3b1SBarry Smith } 68853acd3b1SBarry Smith 68953acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/ 69053acd3b1SBarry Smith #undef __FUNCT__ 691e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private" 69253acd3b1SBarry Smith /*@C 69353acd3b1SBarry Smith PetscOptionsInt - Gets the integer value for a particular option in the database. 69453acd3b1SBarry Smith 6953f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 69653acd3b1SBarry Smith 69753acd3b1SBarry Smith Input Parameters: 69853acd3b1SBarry Smith + opt - option name 69953acd3b1SBarry Smith . text - short string that describes the option 70053acd3b1SBarry Smith . man - manual page with additional information on option 7010fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 7020fdccdaeSBarry Smith $ PetscOptionsInt(..., obj->value,&object->value,...) or 7030fdccdaeSBarry Smith $ value = defaultvalue 7040fdccdaeSBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 7050fdccdaeSBarry Smith $ if (flg) { 70653acd3b1SBarry Smith 70753acd3b1SBarry Smith Output Parameter: 70853acd3b1SBarry Smith + value - the integer value to return 70953acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 71053acd3b1SBarry Smith 71153acd3b1SBarry Smith Level: beginner 71253acd3b1SBarry Smith 71353acd3b1SBarry Smith Concepts: options database^has int 71453acd3b1SBarry Smith 71553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 71653acd3b1SBarry Smith 71753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 718acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 719acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 72053acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 72153acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 722acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 723a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 72453acd3b1SBarry Smith @*/ 7258c34d3f5SBarry Smith PetscErrorCode PetscOptionsInt_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *set) 72653acd3b1SBarry Smith { 72753acd3b1SBarry Smith PetscErrorCode ierr; 7288c34d3f5SBarry Smith PetscOption amsopt; 72912655325SBarry Smith PetscBool wasset; 73053acd3b1SBarry Smith 73153acd3b1SBarry Smith PetscFunctionBegin; 732e55864a3SBarry Smith if (!PetscOptionsObject->count) { 733e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr); 7346356e834SBarry Smith ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr); 73512655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 7363e211508SBarry Smith 73712655325SBarry Smith ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,¤tvalue,&wasset);CHKERRQ(ierr); 7383e211508SBarry Smith if (wasset) { 73912655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 7403e211508SBarry Smith } 741af6d86caSBarry Smith } 742e55864a3SBarry Smith ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 743e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 7441a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 74553acd3b1SBarry Smith } 74653acd3b1SBarry Smith PetscFunctionReturn(0); 74753acd3b1SBarry Smith } 74853acd3b1SBarry Smith 74953acd3b1SBarry Smith #undef __FUNCT__ 750e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private" 75153acd3b1SBarry Smith /*@C 75253acd3b1SBarry Smith PetscOptionsString - Gets the string value for a particular option in the database. 75353acd3b1SBarry Smith 7543f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 75553acd3b1SBarry Smith 75653acd3b1SBarry Smith Input Parameters: 75753acd3b1SBarry Smith + opt - option name 75853acd3b1SBarry Smith . text - short string that describes the option 75953acd3b1SBarry Smith . man - manual page with additional information on option 7600fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 761bcbf2dc5SJed Brown - len - length of the result string including null terminator 76253acd3b1SBarry Smith 76353acd3b1SBarry Smith Output Parameter: 76453acd3b1SBarry Smith + value - the value to return 76553acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 76653acd3b1SBarry Smith 76753acd3b1SBarry Smith Level: beginner 76853acd3b1SBarry Smith 76953acd3b1SBarry Smith Concepts: options database^has int 77053acd3b1SBarry Smith 77153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 77253acd3b1SBarry Smith 7737fccdfe4SBarry 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). 7747fccdfe4SBarry Smith 77553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 776acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 777acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(), 77853acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 77953acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 780acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 781a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 78253acd3b1SBarry Smith @*/ 7838c34d3f5SBarry Smith PetscErrorCode PetscOptionsString_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool *set) 78453acd3b1SBarry Smith { 78553acd3b1SBarry Smith PetscErrorCode ierr; 7868c34d3f5SBarry Smith PetscOption amsopt; 78753acd3b1SBarry Smith 78853acd3b1SBarry Smith PetscFunctionBegin; 7891a1499c8SBarry Smith if (!PetscOptionsObject->count) { 7901a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 79164facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 7920fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 793af6d86caSBarry Smith } 794e55864a3SBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 795e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 7961a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 79753acd3b1SBarry Smith } 79853acd3b1SBarry Smith PetscFunctionReturn(0); 79953acd3b1SBarry Smith } 80053acd3b1SBarry Smith 80153acd3b1SBarry Smith #undef __FUNCT__ 802e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private" 80353acd3b1SBarry Smith /*@C 80453acd3b1SBarry Smith PetscOptionsReal - Gets the PetscReal value for a particular option in the database. 80553acd3b1SBarry Smith 8063f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 80753acd3b1SBarry Smith 80853acd3b1SBarry Smith Input Parameters: 80953acd3b1SBarry Smith + opt - option name 81053acd3b1SBarry Smith . text - short string that describes the option 81153acd3b1SBarry Smith . man - manual page with additional information on option 8120fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 8130fdccdaeSBarry Smith $ PetscOptionsReal(..., obj->value,&object->value,...) or 8140fdccdaeSBarry Smith $ value = defaultvalue 8150fdccdaeSBarry Smith $ PetscOptionsReal(..., value,&value,&flg); 8160fdccdaeSBarry Smith $ if (flg) { 81753acd3b1SBarry Smith 81853acd3b1SBarry Smith Output Parameter: 81953acd3b1SBarry Smith + value - the value to return 82053acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 82153acd3b1SBarry Smith 82253acd3b1SBarry Smith Level: beginner 82353acd3b1SBarry Smith 82453acd3b1SBarry Smith Concepts: options database^has int 82553acd3b1SBarry Smith 82653acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 82753acd3b1SBarry Smith 82853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 829acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 830acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 83153acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 83253acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 833acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 834a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 83553acd3b1SBarry Smith @*/ 8368c34d3f5SBarry Smith PetscErrorCode PetscOptionsReal_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 83753acd3b1SBarry Smith { 83853acd3b1SBarry Smith PetscErrorCode ierr; 8398c34d3f5SBarry Smith PetscOption amsopt; 84053acd3b1SBarry Smith 84153acd3b1SBarry Smith PetscFunctionBegin; 842e55864a3SBarry Smith if (!PetscOptionsObject->count) { 843e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr); 844538aa990SBarry Smith ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 845a297a907SKarl Rupp 8460fdccdaeSBarry Smith *(PetscReal*)amsopt->data = currentvalue; 847538aa990SBarry Smith } 8481a1499c8SBarry Smith ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 8491a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 8501a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr); 85153acd3b1SBarry Smith } 85253acd3b1SBarry Smith PetscFunctionReturn(0); 85353acd3b1SBarry Smith } 85453acd3b1SBarry Smith 85553acd3b1SBarry Smith #undef __FUNCT__ 856e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private" 85753acd3b1SBarry Smith /*@C 85853acd3b1SBarry Smith PetscOptionsScalar - Gets the scalar value for a particular option in the database. 85953acd3b1SBarry Smith 8603f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 86153acd3b1SBarry Smith 86253acd3b1SBarry Smith Input Parameters: 86353acd3b1SBarry Smith + opt - option name 86453acd3b1SBarry Smith . text - short string that describes the option 86553acd3b1SBarry Smith . man - manual page with additional information on option 8660fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 8670fdccdaeSBarry Smith $ PetscOptionsScalar(..., obj->value,&object->value,...) or 8680fdccdaeSBarry Smith $ value = defaultvalue 8690fdccdaeSBarry Smith $ PetscOptionsScalar(..., value,&value,&flg); 8700fdccdaeSBarry Smith $ if (flg) { 8710fdccdaeSBarry Smith 87253acd3b1SBarry Smith 87353acd3b1SBarry Smith Output Parameter: 87453acd3b1SBarry Smith + value - the value to return 87553acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 87653acd3b1SBarry Smith 87753acd3b1SBarry Smith Level: beginner 87853acd3b1SBarry Smith 87953acd3b1SBarry Smith Concepts: options database^has int 88053acd3b1SBarry Smith 88153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 88253acd3b1SBarry Smith 88353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 884acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 885acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 88653acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 88753acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 888acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 889a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 89053acd3b1SBarry Smith @*/ 8918c34d3f5SBarry Smith PetscErrorCode PetscOptionsScalar_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 89253acd3b1SBarry Smith { 89353acd3b1SBarry Smith PetscErrorCode ierr; 89453acd3b1SBarry Smith 89553acd3b1SBarry Smith PetscFunctionBegin; 89653acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 8970fdccdaeSBarry Smith ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr); 89853acd3b1SBarry Smith #else 899e55864a3SBarry Smith ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 90053acd3b1SBarry Smith #endif 90153acd3b1SBarry Smith PetscFunctionReturn(0); 90253acd3b1SBarry Smith } 90353acd3b1SBarry Smith 90453acd3b1SBarry Smith #undef __FUNCT__ 905e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private" 90653acd3b1SBarry Smith /*@C 90790d69ab7SBarry 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 90890d69ab7SBarry Smith its value is set to false. 90953acd3b1SBarry Smith 9103f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 91153acd3b1SBarry Smith 91253acd3b1SBarry Smith Input Parameters: 91353acd3b1SBarry Smith + opt - option name 91453acd3b1SBarry Smith . text - short string that describes the option 91553acd3b1SBarry Smith - man - manual page with additional information on option 91653acd3b1SBarry Smith 91753acd3b1SBarry Smith Output Parameter: 91853acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 91953acd3b1SBarry Smith 92053acd3b1SBarry Smith Level: beginner 92153acd3b1SBarry Smith 92253acd3b1SBarry Smith Concepts: options database^has int 92353acd3b1SBarry Smith 92453acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 92553acd3b1SBarry Smith 92653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 927acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 928acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 92953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 93053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 931acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 932a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 93353acd3b1SBarry Smith @*/ 9348c34d3f5SBarry Smith PetscErrorCode PetscOptionsName_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 93553acd3b1SBarry Smith { 93653acd3b1SBarry Smith PetscErrorCode ierr; 9378c34d3f5SBarry Smith PetscOption amsopt; 93853acd3b1SBarry Smith 93953acd3b1SBarry Smith PetscFunctionBegin; 940e55864a3SBarry Smith if (!PetscOptionsObject->count) { 941e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 942ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 943a297a907SKarl Rupp 944ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 9451ae3d29cSBarry Smith } 946e55864a3SBarry Smith ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr); 947e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 948e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 94953acd3b1SBarry Smith } 95053acd3b1SBarry Smith PetscFunctionReturn(0); 95153acd3b1SBarry Smith } 95253acd3b1SBarry Smith 95353acd3b1SBarry Smith #undef __FUNCT__ 954e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private" 95553acd3b1SBarry Smith /*@C 956a264d7a6SBarry Smith PetscOptionsFList - Puts a list of option values that a single one may be selected from 95753acd3b1SBarry Smith 9583f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 95953acd3b1SBarry Smith 96053acd3b1SBarry Smith Input Parameters: 96153acd3b1SBarry Smith + opt - option name 96253acd3b1SBarry Smith . text - short string that describes the option 96353acd3b1SBarry Smith . man - manual page with additional information on option 96453acd3b1SBarry Smith . list - the possible choices 9650fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 9660fdccdaeSBarry Smith $ PetscOptionsFlist(..., obj->value,value,len,&flg); 9670fdccdaeSBarry Smith $ if (flg) { 9683cc1e11dSBarry Smith - len - the length of the character array value 96953acd3b1SBarry Smith 97053acd3b1SBarry Smith Output Parameter: 97153acd3b1SBarry Smith + value - the value to return 97253acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 97353acd3b1SBarry Smith 97453acd3b1SBarry Smith Level: intermediate 97553acd3b1SBarry Smith 97653acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 97753acd3b1SBarry Smith 97853acd3b1SBarry Smith See PetscOptionsEList() for when the choices are given in a string array 97953acd3b1SBarry Smith 98053acd3b1SBarry Smith To get a listing of all currently specified options, 98188c29154SBarry Smith see PetscOptionsView() or PetscOptionsGetAll() 98253acd3b1SBarry Smith 983eabe10d7SBarry Smith Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list 984eabe10d7SBarry Smith 98553acd3b1SBarry Smith Concepts: options database^list 98653acd3b1SBarry Smith 98753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 988acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 98953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 99053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 991acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 992a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum() 99353acd3b1SBarry Smith @*/ 9948c34d3f5SBarry Smith PetscErrorCode PetscOptionsFList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool *set) 99553acd3b1SBarry Smith { 99653acd3b1SBarry Smith PetscErrorCode ierr; 9978c34d3f5SBarry Smith PetscOption amsopt; 99853acd3b1SBarry Smith 99953acd3b1SBarry Smith PetscFunctionBegin; 10001a1499c8SBarry Smith if (!PetscOptionsObject->count) { 10011a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr); 100264facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 10030fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 10043cc1e11dSBarry Smith amsopt->flist = list; 10053cc1e11dSBarry Smith } 10061a1499c8SBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 10071a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10081a1499c8SBarry Smith ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr); 100953acd3b1SBarry Smith } 101053acd3b1SBarry Smith PetscFunctionReturn(0); 101153acd3b1SBarry Smith } 101253acd3b1SBarry Smith 101353acd3b1SBarry Smith #undef __FUNCT__ 1014e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private" 101553acd3b1SBarry Smith /*@C 101653acd3b1SBarry Smith PetscOptionsEList - Puts a list of option values that a single one may be selected from 101753acd3b1SBarry Smith 10183f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 101953acd3b1SBarry Smith 102053acd3b1SBarry Smith Input Parameters: 102153acd3b1SBarry Smith + opt - option name 102253acd3b1SBarry Smith . ltext - short string that describes the option 102353acd3b1SBarry Smith . man - manual page with additional information on option 1024a264d7a6SBarry Smith . list - the possible choices (one of these must be selected, anything else is invalid) 102553acd3b1SBarry Smith . ntext - number of choices 10260fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10270fdccdaeSBarry Smith $ PetscOptionsElist(..., obj->value,&value,&flg); 10280fdccdaeSBarry Smith $ if (flg) { 10290fdccdaeSBarry Smith 103053acd3b1SBarry Smith 103153acd3b1SBarry Smith Output Parameter: 103253acd3b1SBarry Smith + value - the index of the value to return 103353acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 103453acd3b1SBarry Smith 103553acd3b1SBarry Smith Level: intermediate 103653acd3b1SBarry Smith 103753acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 103853acd3b1SBarry Smith 1039a264d7a6SBarry Smith See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 104053acd3b1SBarry Smith 104153acd3b1SBarry Smith Concepts: options database^list 104253acd3b1SBarry Smith 104353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1044acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 104553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 104653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1047acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1048a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEnum() 104953acd3b1SBarry Smith @*/ 10508c34d3f5SBarry Smith PetscErrorCode PetscOptionsEList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool *set) 105153acd3b1SBarry Smith { 105253acd3b1SBarry Smith PetscErrorCode ierr; 105353acd3b1SBarry Smith PetscInt i; 10548c34d3f5SBarry Smith PetscOption amsopt; 105553acd3b1SBarry Smith 105653acd3b1SBarry Smith PetscFunctionBegin; 10571a1499c8SBarry Smith if (!PetscOptionsObject->count) { 10581a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr); 105964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 10600fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 10611ae3d29cSBarry Smith amsopt->list = list; 10621ae3d29cSBarry Smith amsopt->nlist = ntext; 10631ae3d29cSBarry Smith } 10641a1499c8SBarry Smith ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr); 10651a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10661a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr); 106753acd3b1SBarry Smith for (i=0; i<ntext; i++) { 1068e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr); 106953acd3b1SBarry Smith } 1070e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr); 107153acd3b1SBarry Smith } 107253acd3b1SBarry Smith PetscFunctionReturn(0); 107353acd3b1SBarry Smith } 107453acd3b1SBarry Smith 107553acd3b1SBarry Smith #undef __FUNCT__ 1076e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private" 107753acd3b1SBarry Smith /*@C 1078acfcf0e5SJed Brown PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 1079d5649816SBarry Smith which at most a single value can be true. 108053acd3b1SBarry Smith 10813f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 108253acd3b1SBarry Smith 108353acd3b1SBarry Smith Input Parameters: 108453acd3b1SBarry Smith + opt - option name 108553acd3b1SBarry Smith . text - short string that describes the option 108653acd3b1SBarry Smith - man - manual page with additional information on option 108753acd3b1SBarry Smith 108853acd3b1SBarry Smith Output Parameter: 108953acd3b1SBarry Smith . flg - whether that option was set or not 109053acd3b1SBarry Smith 109153acd3b1SBarry Smith Level: intermediate 109253acd3b1SBarry Smith 109353acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 109453acd3b1SBarry Smith 1095acfcf0e5SJed Brown Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd() 109653acd3b1SBarry Smith 109753acd3b1SBarry Smith Concepts: options database^logical group 109853acd3b1SBarry Smith 109953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1100acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 110153acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 110253acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1103acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1104a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 110553acd3b1SBarry Smith @*/ 11068c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 110753acd3b1SBarry Smith { 110853acd3b1SBarry Smith PetscErrorCode ierr; 11098c34d3f5SBarry Smith PetscOption amsopt; 111053acd3b1SBarry Smith 111153acd3b1SBarry Smith PetscFunctionBegin; 1112e55864a3SBarry Smith if (!PetscOptionsObject->count) { 111383355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1114ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1115a297a907SKarl Rupp 1116ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 11171ae3d29cSBarry Smith } 111868b16fdaSBarry Smith *flg = PETSC_FALSE; 1119e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1120e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1121e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," Pick at most one of -------------\n");CHKERRQ(ierr); 1122e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 112353acd3b1SBarry Smith } 112453acd3b1SBarry Smith PetscFunctionReturn(0); 112553acd3b1SBarry Smith } 112653acd3b1SBarry Smith 112753acd3b1SBarry Smith #undef __FUNCT__ 1128e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private" 112953acd3b1SBarry Smith /*@C 1130acfcf0e5SJed Brown PetscOptionsBoolGroup - One in a series of logical queries on the options database for 1131d5649816SBarry Smith which at most a single value can be true. 113253acd3b1SBarry Smith 11333f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 113453acd3b1SBarry Smith 113553acd3b1SBarry Smith Input Parameters: 113653acd3b1SBarry Smith + opt - option name 113753acd3b1SBarry Smith . text - short string that describes the option 113853acd3b1SBarry Smith - man - manual page with additional information on option 113953acd3b1SBarry Smith 114053acd3b1SBarry Smith Output Parameter: 114153acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 114253acd3b1SBarry Smith 114353acd3b1SBarry Smith Level: intermediate 114453acd3b1SBarry Smith 114553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 114653acd3b1SBarry Smith 1147acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd() 114853acd3b1SBarry Smith 114953acd3b1SBarry Smith Concepts: options database^logical group 115053acd3b1SBarry Smith 115153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1152acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 115353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 115453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1155acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1156a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 115753acd3b1SBarry Smith @*/ 11588c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 115953acd3b1SBarry Smith { 116053acd3b1SBarry Smith PetscErrorCode ierr; 11618c34d3f5SBarry Smith PetscOption amsopt; 116253acd3b1SBarry Smith 116353acd3b1SBarry Smith PetscFunctionBegin; 1164e55864a3SBarry Smith if (!PetscOptionsObject->count) { 116583355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1166ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1167a297a907SKarl Rupp 1168ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 11691ae3d29cSBarry Smith } 117017326d04SJed Brown *flg = PETSC_FALSE; 1171e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1172e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1173e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 117453acd3b1SBarry Smith } 117553acd3b1SBarry Smith PetscFunctionReturn(0); 117653acd3b1SBarry Smith } 117753acd3b1SBarry Smith 117853acd3b1SBarry Smith #undef __FUNCT__ 1179e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private" 118053acd3b1SBarry Smith /*@C 1181acfcf0e5SJed Brown PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 1182d5649816SBarry Smith which at most a single value can be true. 118353acd3b1SBarry Smith 11843f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 118553acd3b1SBarry Smith 118653acd3b1SBarry Smith Input Parameters: 118753acd3b1SBarry Smith + opt - option name 118853acd3b1SBarry Smith . text - short string that describes the option 118953acd3b1SBarry Smith - man - manual page with additional information on option 119053acd3b1SBarry Smith 119153acd3b1SBarry Smith Output Parameter: 119253acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 119353acd3b1SBarry Smith 119453acd3b1SBarry Smith Level: intermediate 119553acd3b1SBarry Smith 119653acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 119753acd3b1SBarry Smith 1198acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() 119953acd3b1SBarry Smith 120053acd3b1SBarry Smith Concepts: options database^logical group 120153acd3b1SBarry Smith 120253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1203acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 120453acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 120553acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1206acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1207a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 120853acd3b1SBarry Smith @*/ 12098c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 121053acd3b1SBarry Smith { 121153acd3b1SBarry Smith PetscErrorCode ierr; 12128c34d3f5SBarry Smith PetscOption amsopt; 121353acd3b1SBarry Smith 121453acd3b1SBarry Smith PetscFunctionBegin; 1215e55864a3SBarry Smith if (!PetscOptionsObject->count) { 121683355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1217ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1218a297a907SKarl Rupp 1219ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 12201ae3d29cSBarry Smith } 122117326d04SJed Brown *flg = PETSC_FALSE; 1222e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1223e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1224e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 122553acd3b1SBarry Smith } 122653acd3b1SBarry Smith PetscFunctionReturn(0); 122753acd3b1SBarry Smith } 122853acd3b1SBarry Smith 122953acd3b1SBarry Smith #undef __FUNCT__ 1230e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private" 123153acd3b1SBarry Smith /*@C 1232acfcf0e5SJed Brown PetscOptionsBool - Determines if a particular option is in the database with a true or false 123353acd3b1SBarry Smith 12343f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 123553acd3b1SBarry Smith 123653acd3b1SBarry Smith Input Parameters: 123753acd3b1SBarry Smith + opt - option name 123853acd3b1SBarry Smith . text - short string that describes the option 1239868c398cSBarry Smith . man - manual page with additional information on option 124094ae4db5SBarry Smith - currentvalue - the current value 124153acd3b1SBarry Smith 124253acd3b1SBarry Smith Output Parameter: 124353acd3b1SBarry Smith . flg - PETSC_TRUE or PETSC_FALSE 124453acd3b1SBarry Smith . set - PETSC_TRUE if found, else PETSC_FALSE 124553acd3b1SBarry Smith 124653acd3b1SBarry Smith Level: beginner 124753acd3b1SBarry Smith 124853acd3b1SBarry Smith Concepts: options database^logical 124953acd3b1SBarry Smith 125053acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 125153acd3b1SBarry Smith 125253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 1253acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1254acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 125553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 125653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1257acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1258a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 125953acd3b1SBarry Smith @*/ 12608c34d3f5SBarry Smith PetscErrorCode PetscOptionsBool_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 126153acd3b1SBarry Smith { 126253acd3b1SBarry Smith PetscErrorCode ierr; 1263ace3abfcSBarry Smith PetscBool iset; 12648c34d3f5SBarry Smith PetscOption amsopt; 126553acd3b1SBarry Smith 126653acd3b1SBarry Smith PetscFunctionBegin; 1267e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1268e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1269ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1270a297a907SKarl Rupp 127194ae4db5SBarry Smith *(PetscBool*)amsopt->data = currentvalue; 1272af6d86caSBarry Smith } 12731a1499c8SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr); 127453acd3b1SBarry Smith if (set) *set = iset; 12751a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 127694ae4db5SBarry Smith const char *v = PetscBools[currentvalue]; 12771a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr); 127853acd3b1SBarry Smith } 127953acd3b1SBarry Smith PetscFunctionReturn(0); 128053acd3b1SBarry Smith } 128153acd3b1SBarry Smith 128253acd3b1SBarry Smith #undef __FUNCT__ 1283e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private" 128453acd3b1SBarry Smith /*@C 128553acd3b1SBarry Smith PetscOptionsRealArray - Gets an array of double values for a particular 128653acd3b1SBarry Smith option in the database. The values must be separated with commas with 128753acd3b1SBarry Smith no intervening spaces. 128853acd3b1SBarry Smith 12893f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 129053acd3b1SBarry Smith 129153acd3b1SBarry Smith Input Parameters: 129253acd3b1SBarry Smith + opt - the option one is seeking 129353acd3b1SBarry Smith . text - short string describing option 129453acd3b1SBarry Smith . man - manual page for option 129553acd3b1SBarry Smith - nmax - maximum number of values 129653acd3b1SBarry Smith 129753acd3b1SBarry Smith Output Parameter: 129853acd3b1SBarry Smith + value - location to copy values 129953acd3b1SBarry Smith . nmax - actual number of values found 130053acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 130153acd3b1SBarry Smith 130253acd3b1SBarry Smith Level: beginner 130353acd3b1SBarry Smith 130453acd3b1SBarry Smith Notes: 130553acd3b1SBarry Smith The user should pass in an array of doubles 130653acd3b1SBarry Smith 130753acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 130853acd3b1SBarry Smith 130953acd3b1SBarry Smith Concepts: options database^array of strings 131053acd3b1SBarry Smith 131153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1312acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 131353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 131453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1315acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1316a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 131753acd3b1SBarry Smith @*/ 13188c34d3f5SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 131953acd3b1SBarry Smith { 132053acd3b1SBarry Smith PetscErrorCode ierr; 132153acd3b1SBarry Smith PetscInt i; 13228c34d3f5SBarry Smith PetscOption amsopt; 132353acd3b1SBarry Smith 132453acd3b1SBarry Smith PetscFunctionBegin; 1325e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1326e26ddf31SBarry Smith PetscReal *vals; 1327e26ddf31SBarry Smith 1328e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr); 1329e55864a3SBarry Smith ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 1330e26ddf31SBarry Smith vals = (PetscReal*)amsopt->data; 1331e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1332e26ddf31SBarry Smith amsopt->arraylength = *n; 1333e26ddf31SBarry Smith } 1334e55864a3SBarry Smith ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1335e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1336a519f713SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr); 133753acd3b1SBarry Smith for (i=1; i<*n; i++) { 1338e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr); 133953acd3b1SBarry Smith } 1340e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 134153acd3b1SBarry Smith } 134253acd3b1SBarry Smith PetscFunctionReturn(0); 134353acd3b1SBarry Smith } 134453acd3b1SBarry Smith 134553acd3b1SBarry Smith 134653acd3b1SBarry Smith #undef __FUNCT__ 1347e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private" 134853acd3b1SBarry Smith /*@C 134953acd3b1SBarry Smith PetscOptionsIntArray - Gets an array of integers for a particular 1350b32a342fSShri Abhyankar option in the database. 135153acd3b1SBarry Smith 13523f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 135353acd3b1SBarry Smith 135453acd3b1SBarry Smith Input Parameters: 135553acd3b1SBarry Smith + opt - the option one is seeking 135653acd3b1SBarry Smith . text - short string describing option 135753acd3b1SBarry Smith . man - manual page for option 1358f8a50e2bSBarry Smith - n - maximum number of values 135953acd3b1SBarry Smith 136053acd3b1SBarry Smith Output Parameter: 136153acd3b1SBarry Smith + value - location to copy values 1362f8a50e2bSBarry Smith . n - actual number of values found 136353acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 136453acd3b1SBarry Smith 136553acd3b1SBarry Smith Level: beginner 136653acd3b1SBarry Smith 136753acd3b1SBarry Smith Notes: 1368b32a342fSShri Abhyankar The array can be passed as 1369b32a342fSShri Abhyankar a comma seperated list: 0,1,2,3,4,5,6,7 13700fd488f5SShri Abhyankar a range (start-end+1): 0-8 13710fd488f5SShri Abhyankar a range with given increment (start-end+1:inc): 0-7:2 13720fd488f5SShri Abhyankar a combination of values and ranges seperated by commas: 0,1-8,8-15:2 1373b32a342fSShri Abhyankar 1374b32a342fSShri Abhyankar There must be no intervening spaces between the values. 137553acd3b1SBarry Smith 137653acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 137753acd3b1SBarry Smith 1378b32a342fSShri Abhyankar Concepts: options database^array of ints 137953acd3b1SBarry Smith 138053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1381acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 138253acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 138353acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1384acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1385a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray() 138653acd3b1SBarry Smith @*/ 13878c34d3f5SBarry Smith PetscErrorCode PetscOptionsIntArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 138853acd3b1SBarry Smith { 138953acd3b1SBarry Smith PetscErrorCode ierr; 139053acd3b1SBarry Smith PetscInt i; 13918c34d3f5SBarry Smith PetscOption amsopt; 139253acd3b1SBarry Smith 139353acd3b1SBarry Smith PetscFunctionBegin; 1394e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1395e26ddf31SBarry Smith PetscInt *vals; 1396e26ddf31SBarry Smith 1397e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr); 1398854ce69bSBarry Smith ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr); 1399e26ddf31SBarry Smith vals = (PetscInt*)amsopt->data; 1400e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1401e26ddf31SBarry Smith amsopt->arraylength = *n; 1402e26ddf31SBarry Smith } 1403e55864a3SBarry Smith ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1404e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1405e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 140653acd3b1SBarry Smith for (i=1; i<*n; i++) { 1407e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 140853acd3b1SBarry Smith } 1409e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 141053acd3b1SBarry Smith } 141153acd3b1SBarry Smith PetscFunctionReturn(0); 141253acd3b1SBarry Smith } 141353acd3b1SBarry Smith 141453acd3b1SBarry Smith #undef __FUNCT__ 1415e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private" 141653acd3b1SBarry Smith /*@C 141753acd3b1SBarry Smith PetscOptionsStringArray - Gets an array of string values for a particular 141853acd3b1SBarry Smith option in the database. The values must be separated with commas with 141953acd3b1SBarry Smith no intervening spaces. 142053acd3b1SBarry Smith 14213f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 142253acd3b1SBarry Smith 142353acd3b1SBarry Smith Input Parameters: 142453acd3b1SBarry Smith + opt - the option one is seeking 142553acd3b1SBarry Smith . text - short string describing option 142653acd3b1SBarry Smith . man - manual page for option 142753acd3b1SBarry Smith - nmax - maximum number of strings 142853acd3b1SBarry Smith 142953acd3b1SBarry Smith Output Parameter: 143053acd3b1SBarry Smith + value - location to copy strings 143153acd3b1SBarry Smith . nmax - actual number of strings found 143253acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 143353acd3b1SBarry Smith 143453acd3b1SBarry Smith Level: beginner 143553acd3b1SBarry Smith 143653acd3b1SBarry Smith Notes: 143753acd3b1SBarry Smith The user should pass in an array of pointers to char, to hold all the 143853acd3b1SBarry Smith strings returned by this function. 143953acd3b1SBarry Smith 144053acd3b1SBarry Smith The user is responsible for deallocating the strings that are 144153acd3b1SBarry Smith returned. The Fortran interface for this routine is not supported. 144253acd3b1SBarry Smith 144353acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 144453acd3b1SBarry Smith 144553acd3b1SBarry Smith Concepts: options database^array of strings 144653acd3b1SBarry Smith 144753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1448acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 144953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 145053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1451acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1452a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 145353acd3b1SBarry Smith @*/ 14548c34d3f5SBarry Smith PetscErrorCode PetscOptionsStringArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 145553acd3b1SBarry Smith { 145653acd3b1SBarry Smith PetscErrorCode ierr; 14578c34d3f5SBarry Smith PetscOption amsopt; 145853acd3b1SBarry Smith 145953acd3b1SBarry Smith PetscFunctionBegin; 1460e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1461e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr); 1462854ce69bSBarry Smith ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr); 1463a297a907SKarl Rupp 14641ae3d29cSBarry Smith amsopt->arraylength = *nmax; 14651ae3d29cSBarry Smith } 1466e55864a3SBarry Smith ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr); 1467e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1468e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 146953acd3b1SBarry Smith } 147053acd3b1SBarry Smith PetscFunctionReturn(0); 147153acd3b1SBarry Smith } 147253acd3b1SBarry Smith 1473e2446a98SMatthew Knepley #undef __FUNCT__ 1474e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private" 1475e2446a98SMatthew Knepley /*@C 1476acfcf0e5SJed Brown PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 1477e2446a98SMatthew Knepley option in the database. The values must be separated with commas with 1478e2446a98SMatthew Knepley no intervening spaces. 1479e2446a98SMatthew Knepley 14803f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 1481e2446a98SMatthew Knepley 1482e2446a98SMatthew Knepley Input Parameters: 1483e2446a98SMatthew Knepley + opt - the option one is seeking 1484e2446a98SMatthew Knepley . text - short string describing option 1485e2446a98SMatthew Knepley . man - manual page for option 1486e2446a98SMatthew Knepley - nmax - maximum number of values 1487e2446a98SMatthew Knepley 1488e2446a98SMatthew Knepley Output Parameter: 1489e2446a98SMatthew Knepley + value - location to copy values 1490e2446a98SMatthew Knepley . nmax - actual number of values found 1491e2446a98SMatthew Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 1492e2446a98SMatthew Knepley 1493e2446a98SMatthew Knepley Level: beginner 1494e2446a98SMatthew Knepley 1495e2446a98SMatthew Knepley Notes: 1496e2446a98SMatthew Knepley The user should pass in an array of doubles 1497e2446a98SMatthew Knepley 1498e2446a98SMatthew Knepley Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1499e2446a98SMatthew Knepley 1500e2446a98SMatthew Knepley Concepts: options database^array of strings 1501e2446a98SMatthew Knepley 1502e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1503acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1504e2446a98SMatthew Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1505e2446a98SMatthew Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1506acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1507a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 1508e2446a98SMatthew Knepley @*/ 15098c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 1510e2446a98SMatthew Knepley { 1511e2446a98SMatthew Knepley PetscErrorCode ierr; 1512e2446a98SMatthew Knepley PetscInt i; 15138c34d3f5SBarry Smith PetscOption amsopt; 1514e2446a98SMatthew Knepley 1515e2446a98SMatthew Knepley PetscFunctionBegin; 1516e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1517ace3abfcSBarry Smith PetscBool *vals; 15181ae3d29cSBarry Smith 151983355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr); 15201a1499c8SBarry Smith ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr); 1521ace3abfcSBarry Smith vals = (PetscBool*)amsopt->data; 15221ae3d29cSBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 15231ae3d29cSBarry Smith amsopt->arraylength = *n; 15241ae3d29cSBarry Smith } 1525e55864a3SBarry Smith ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1526e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1527e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 1528e2446a98SMatthew Knepley for (i=1; i<*n; i++) { 1529e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 1530e2446a98SMatthew Knepley } 1531e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 1532e2446a98SMatthew Knepley } 1533e2446a98SMatthew Knepley PetscFunctionReturn(0); 1534e2446a98SMatthew Knepley } 1535e2446a98SMatthew Knepley 15368cc676e6SMatthew G Knepley #undef __FUNCT__ 1537e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private" 15388cc676e6SMatthew G Knepley /*@C 1539d1da0b69SBarry Smith PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user 15408cc676e6SMatthew G Knepley 15418cc676e6SMatthew G Knepley Logically Collective on the communicator passed in PetscOptionsBegin() 15428cc676e6SMatthew G Knepley 15438cc676e6SMatthew G Knepley Input Parameters: 15448cc676e6SMatthew G Knepley + opt - option name 15458cc676e6SMatthew G Knepley . text - short string that describes the option 15468cc676e6SMatthew G Knepley - man - manual page with additional information on option 15478cc676e6SMatthew G Knepley 15488cc676e6SMatthew G Knepley Output Parameter: 15498cc676e6SMatthew G Knepley + viewer - the viewer 15508cc676e6SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 15518cc676e6SMatthew G Knepley 15528cc676e6SMatthew G Knepley Level: beginner 15538cc676e6SMatthew G Knepley 15548cc676e6SMatthew G Knepley Concepts: options database^has int 15558cc676e6SMatthew G Knepley 15568cc676e6SMatthew G Knepley Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 15578cc676e6SMatthew G Knepley 15585a7113b9SPatrick Sanan See PetscOptionsGetViewer() for the format of the supplied viewer and its options 15598cc676e6SMatthew G Knepley 15608cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 15618cc676e6SMatthew G Knepley PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 15628cc676e6SMatthew G Knepley PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 15638cc676e6SMatthew G Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 15648cc676e6SMatthew G Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 15658cc676e6SMatthew G Knepley PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1566a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 15678cc676e6SMatthew G Knepley @*/ 15688c34d3f5SBarry Smith PetscErrorCode PetscOptionsViewer_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 15698cc676e6SMatthew G Knepley { 15708cc676e6SMatthew G Knepley PetscErrorCode ierr; 15718c34d3f5SBarry Smith PetscOption amsopt; 15728cc676e6SMatthew G Knepley 15738cc676e6SMatthew G Knepley PetscFunctionBegin; 15741a1499c8SBarry Smith if (!PetscOptionsObject->count) { 15751a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 157664facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 15775b02f95dSBarry Smith ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr); 15788cc676e6SMatthew G Knepley } 1579e55864a3SBarry Smith ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr); 1580e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1581e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr); 15828cc676e6SMatthew G Knepley } 15838cc676e6SMatthew G Knepley PetscFunctionReturn(0); 15848cc676e6SMatthew G Knepley } 15858cc676e6SMatthew G Knepley 158653acd3b1SBarry Smith 158753acd3b1SBarry Smith #undef __FUNCT__ 158853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead" 158953acd3b1SBarry Smith /*@C 1590b52f573bSBarry Smith PetscOptionsHead - Puts a heading before listing any more published options. Used, for example, 159153acd3b1SBarry Smith in KSPSetFromOptions_GMRES(). 159253acd3b1SBarry Smith 15933f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 159453acd3b1SBarry Smith 159553acd3b1SBarry Smith Input Parameter: 159653acd3b1SBarry Smith . head - the heading text 159753acd3b1SBarry Smith 159853acd3b1SBarry Smith 159953acd3b1SBarry Smith Level: intermediate 160053acd3b1SBarry Smith 160153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 160253acd3b1SBarry Smith 1603b52f573bSBarry Smith Can be followed by a call to PetscOptionsTail() in the same function. 160453acd3b1SBarry Smith 160553acd3b1SBarry Smith Concepts: options database^subheading 160653acd3b1SBarry Smith 160753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1608acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 160953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 161053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1611acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1612a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 161353acd3b1SBarry Smith @*/ 16148c34d3f5SBarry Smith PetscErrorCode PetscOptionsHead(PetscOptions *PetscOptionsObject,const char head[]) 161553acd3b1SBarry Smith { 161653acd3b1SBarry Smith PetscErrorCode ierr; 161753acd3b1SBarry Smith 161853acd3b1SBarry Smith PetscFunctionBegin; 1619e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1620e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s\n",head);CHKERRQ(ierr); 162153acd3b1SBarry Smith } 162253acd3b1SBarry Smith PetscFunctionReturn(0); 162353acd3b1SBarry Smith } 162453acd3b1SBarry Smith 162553acd3b1SBarry Smith 162653acd3b1SBarry Smith 162753acd3b1SBarry Smith 162853acd3b1SBarry Smith 162953acd3b1SBarry Smith 1630