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 9af0996ceSBarry Smith #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 1992a409bb0SBarry Smith PetscFunctionBegin; 200e55864a3SBarry Smith ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr); 2016356e834SBarry Smith while (next) { 2026356e834SBarry Smith switch (next->type) { 2036356e834SBarry Smith case OPTION_HEAD: 2046356e834SBarry Smith break; 205e26ddf31SBarry Smith case OPTION_INT_ARRAY: 206e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 207e26ddf31SBarry Smith vald = (PetscInt*) next->data; 208e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 209e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr); 210e26ddf31SBarry Smith if (i < next->arraylength-1) { 211e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 212e26ddf31SBarry Smith } 213e26ddf31SBarry Smith } 214e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 215e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 216e26ddf31SBarry Smith if (str[0]) { 217e26ddf31SBarry Smith PetscToken token; 218e26ddf31SBarry Smith PetscInt n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end; 219e26ddf31SBarry Smith size_t len; 220e26ddf31SBarry Smith char *value; 221ace3abfcSBarry Smith PetscBool foundrange; 222e26ddf31SBarry Smith 223e26ddf31SBarry Smith next->set = PETSC_TRUE; 224e26ddf31SBarry Smith value = str; 225e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 226e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 227e26ddf31SBarry Smith while (n < nmax) { 228e26ddf31SBarry Smith if (!value) break; 229e26ddf31SBarry Smith 230e26ddf31SBarry Smith /* look for form d-D where d and D are integers */ 231e26ddf31SBarry Smith foundrange = PETSC_FALSE; 232e26ddf31SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 233e26ddf31SBarry Smith if (value[0] == '-') i=2; 234e26ddf31SBarry Smith else i=1; 235330cf3c9SBarry Smith for (;i<len; i++) { 236e26ddf31SBarry Smith if (value[i] == '-') { 237e32f2f54SBarry Smith if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 238e26ddf31SBarry Smith value[i] = 0; 239cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 240cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 241e32f2f54SBarry 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); 242e32f2f54SBarry 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); 243e26ddf31SBarry Smith for (; start<end; start++) { 244e26ddf31SBarry Smith *dvalue = start; dvalue++;n++; 245e26ddf31SBarry Smith } 246e26ddf31SBarry Smith foundrange = PETSC_TRUE; 247e26ddf31SBarry Smith break; 248e26ddf31SBarry Smith } 249e26ddf31SBarry Smith } 250e26ddf31SBarry Smith if (!foundrange) { 251cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr); 252e26ddf31SBarry Smith dvalue++; 253e26ddf31SBarry Smith n++; 254e26ddf31SBarry Smith } 255e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 256e26ddf31SBarry Smith } 2578c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 258e26ddf31SBarry Smith } 259e26ddf31SBarry Smith break; 260e26ddf31SBarry Smith case OPTION_REAL_ARRAY: 261e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 262e26ddf31SBarry Smith valr = (PetscReal*) next->data; 263e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 264e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr); 265e26ddf31SBarry Smith if (i < next->arraylength-1) { 266e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 267e26ddf31SBarry Smith } 268e26ddf31SBarry Smith } 269e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 270e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 271e26ddf31SBarry Smith if (str[0]) { 272e26ddf31SBarry Smith PetscToken token; 273e26ddf31SBarry Smith PetscInt n = 0,nmax = next->arraylength; 274e26ddf31SBarry Smith PetscReal *dvalue = (PetscReal*)next->data; 275e26ddf31SBarry Smith char *value; 276e26ddf31SBarry Smith 277e26ddf31SBarry Smith next->set = PETSC_TRUE; 278e26ddf31SBarry Smith value = str; 279e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 280e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 281e26ddf31SBarry Smith while (n < nmax) { 282e26ddf31SBarry Smith if (!value) break; 283cfbddea1SSatish Balay ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 284e26ddf31SBarry Smith dvalue++; 285e26ddf31SBarry Smith n++; 286e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 287e26ddf31SBarry Smith } 2888c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 289e26ddf31SBarry Smith } 290e26ddf31SBarry Smith break; 2916356e834SBarry Smith case OPTION_INT: 292e55864a3SBarry 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); 2933fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 2943fc1eb6aSBarry Smith if (str[0]) { 295d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 296d25d7f95SJed Brown long long lid; 297d25d7f95SJed Brown sscanf(str,"%lld",&lid); 2981a1499c8SBarry 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); 299c272547aSJed Brown #else 300d25d7f95SJed Brown long lid; 301d25d7f95SJed Brown sscanf(str,"%ld",&lid); 3021a1499c8SBarry 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); 303c272547aSJed Brown #endif 304a297a907SKarl Rupp 305d25d7f95SJed Brown next->set = PETSC_TRUE; 306d25d7f95SJed Brown *((PetscInt*)next->data) = (PetscInt)lid; 307aee2cecaSBarry Smith } 308aee2cecaSBarry Smith break; 309aee2cecaSBarry Smith case OPTION_REAL: 310e55864a3SBarry 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); 3113fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3123fc1eb6aSBarry Smith if (str[0]) { 313ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 314a4404d99SBarry Smith sscanf(str,"%e",&ir); 315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE) 316aee2cecaSBarry Smith sscanf(str,"%le",&ir); 317ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 318d9822059SBarry Smith ir = strtoflt128(str,0); 319d9822059SBarry Smith #else 320513dbe71SLisandro Dalcin SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type"); 321a4404d99SBarry Smith #endif 322aee2cecaSBarry Smith next->set = PETSC_TRUE; 323aee2cecaSBarry Smith *((PetscReal*)next->data) = ir; 324aee2cecaSBarry Smith } 325aee2cecaSBarry Smith break; 3267781c08eSBarry Smith case OPTION_BOOL: 32783355fc5SBarry 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); 3287781c08eSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3297781c08eSBarry Smith if (str[0]) { 3307781c08eSBarry Smith ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr); 3317781c08eSBarry Smith next->set = PETSC_TRUE; 3327781c08eSBarry Smith *((PetscBool*)next->data) = bid; 3337781c08eSBarry Smith } 3347781c08eSBarry Smith break; 335aee2cecaSBarry Smith case OPTION_STRING: 336e55864a3SBarry 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); 3373fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3383fc1eb6aSBarry Smith if (str[0]) { 339aee2cecaSBarry Smith next->set = PETSC_TRUE; 34064facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3415b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3426356e834SBarry Smith } 3436356e834SBarry Smith break; 344a264d7a6SBarry Smith case OPTION_FLIST: 345e55864a3SBarry Smith ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr); 3463cc1e11dSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3473cc1e11dSBarry Smith if (str[0]) { 348e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_TRUE; 3493cc1e11dSBarry Smith next->set = PETSC_TRUE; 35064facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3515b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3523cc1e11dSBarry Smith } 3533cc1e11dSBarry Smith break; 354b432afdaSMatthew Knepley default: 355b432afdaSMatthew Knepley break; 3566356e834SBarry Smith } 3576356e834SBarry Smith next = next->next; 3586356e834SBarry Smith } 3596356e834SBarry Smith PetscFunctionReturn(0); 3606356e834SBarry Smith } 3616356e834SBarry Smith 362e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 363e04113cfSBarry Smith #include <petscviewersaws.h> 364d5649816SBarry Smith 365d5649816SBarry Smith static int count = 0; 366d5649816SBarry Smith 367b3506946SBarry Smith #undef __FUNCT__ 368e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy" 369e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void) 370d5649816SBarry Smith { 3712657e9d9SBarry Smith PetscFunctionBegin; 372d5649816SBarry Smith PetscFunctionReturn(0); 373d5649816SBarry Smith } 374d5649816SBarry Smith 3759c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n" 37623a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n" 37723a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n" 378d1fc0251SBarry Smith "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n" 37964bbc9efSBarry Smith "<script>\n" 38064bbc9efSBarry Smith "jQuery(document).ready(function() {\n" 38164bbc9efSBarry Smith "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n" 38264bbc9efSBarry Smith "})\n" 38364bbc9efSBarry Smith "</script>\n" 38464bbc9efSBarry Smith "</head>\n"; 3851423471aSBarry Smith 3861423471aSBarry Smith /* Determines the size and style of the scroll region where PETSc options selectable from users are displayed */ 3871423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>"; 38864bbc9efSBarry Smith 389d5649816SBarry Smith #undef __FUNCT__ 3907781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput" 391b3506946SBarry Smith /* 3927781c08eSBarry Smith PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs 393b3506946SBarry Smith 394b3506946SBarry Smith Bugs: 395b3506946SBarry Smith + All processes must traverse through the exact same set of option queries do to the call to PetscScanString() 396b3506946SBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 397b3506946SBarry Smith - Only works for PetscInt == int, PetscReal == double etc 398b3506946SBarry Smith 399b3506946SBarry Smith 400b3506946SBarry Smith */ 401f7b25cf6SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptions *PetscOptionsObject) 402b3506946SBarry Smith { 403b3506946SBarry Smith PetscErrorCode ierr; 4048c34d3f5SBarry Smith PetscOption next = PetscOptionsObject->next; 405d5649816SBarry Smith static int mancount = 0; 406b3506946SBarry Smith char options[16]; 4077aab2a10SBarry Smith PetscBool changedmethod = PETSC_FALSE; 408a23eb982SSurtai Han PetscBool stopasking = PETSC_FALSE; 40988a9752cSBarry Smith char manname[16],textname[16]; 4102657e9d9SBarry Smith char dir[1024]; 411b3506946SBarry Smith 4122a409bb0SBarry Smith PetscFunctionBegin; 413b3506946SBarry Smith /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */ 414b3506946SBarry Smith sprintf(options,"Options_%d",count++); 415a297a907SKarl Rupp 4167325c4a4SBarry Smith PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */ 4171bc75a8dSBarry Smith 418d50831c4SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr); 41983355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING)); 4207781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr); 42183355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING)); 4222657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN)); 423a23eb982SSurtai Han PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN)); 424b3506946SBarry Smith 425b3506946SBarry Smith while (next) { 426d50831c4SBarry Smith sprintf(manname,"_man_%d",mancount); 4272657e9d9SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr); 4287781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING)); 429d50831c4SBarry Smith sprintf(textname,"_text_%d",mancount++); 4307781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr); 4317781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING)); 4329f32e415SBarry Smith 433b3506946SBarry Smith switch (next->type) { 434b3506946SBarry Smith case OPTION_HEAD: 435b3506946SBarry Smith break; 436b3506946SBarry Smith case OPTION_INT_ARRAY: 4377781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4382657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT)); 439b3506946SBarry Smith break; 440b3506946SBarry Smith case OPTION_REAL_ARRAY: 4417781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4422657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE)); 443b3506946SBarry Smith break; 444b3506946SBarry Smith case OPTION_INT: 4457781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4462657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT)); 447b3506946SBarry Smith break; 448b3506946SBarry Smith case OPTION_REAL: 4497781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4502657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE)); 451b3506946SBarry Smith break; 4527781c08eSBarry Smith case OPTION_BOOL: 4537781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4542657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN)); 4551ae3d29cSBarry Smith break; 4567781c08eSBarry Smith case OPTION_BOOL_ARRAY: 4577781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4582657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN)); 45971f08665SBarry Smith break; 460b3506946SBarry Smith case OPTION_STRING: 4617781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4627781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4631ae3d29cSBarry Smith break; 4641ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 4657781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4662657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING)); 467b3506946SBarry Smith break; 468a264d7a6SBarry Smith case OPTION_FLIST: 469a264d7a6SBarry Smith { 470a264d7a6SBarry Smith PetscInt ntext; 4717781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4727781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 473a264d7a6SBarry Smith ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr); 474a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 475a264d7a6SBarry Smith } 476a264d7a6SBarry Smith break; 4771ae3d29cSBarry Smith case OPTION_ELIST: 478a264d7a6SBarry Smith { 479a264d7a6SBarry Smith PetscInt ntext = next->nlist; 4807781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4817781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 482ead66b60SBarry Smith ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr); 4831ae3d29cSBarry Smith ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr); 484a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 485a264d7a6SBarry Smith } 486a264d7a6SBarry Smith break; 487b3506946SBarry Smith default: 488b3506946SBarry Smith break; 489b3506946SBarry Smith } 490b3506946SBarry Smith next = next->next; 491b3506946SBarry Smith } 492b3506946SBarry Smith 493b3506946SBarry Smith /* wait until accessor has unlocked the memory */ 49464bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader)); 49564bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom)); 4967aab2a10SBarry Smith ierr = PetscSAWsBlock();CHKERRQ(ierr); 49764bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Header,("index.html")); 49864bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2)); 499b3506946SBarry Smith 50088a9752cSBarry Smith /* determine if any values have been set in GUI */ 50183355fc5SBarry Smith next = PetscOptionsObject->next; 50288a9752cSBarry Smith while (next) { 50388a9752cSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 504f7b25cf6SBarry Smith PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set)); 50588a9752cSBarry Smith next = next->next; 50688a9752cSBarry Smith } 50788a9752cSBarry Smith 508b3506946SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 509f7b25cf6SBarry Smith if (changedmethod) PetscOptionsObject->count = -2; 510b3506946SBarry Smith 511a23eb982SSurtai Han if (stopasking) { 512a23eb982SSurtai Han PetscOptionsPublish = PETSC_FALSE; 51312655325SBarry Smith PetscOptionsObject->count = 0;//do not ask for same thing again 514a23eb982SSurtai Han } 515a23eb982SSurtai Han 5169a492a5cSBarry Smith PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options")); 517b3506946SBarry Smith PetscFunctionReturn(0); 518b3506946SBarry Smith } 519b3506946SBarry Smith #endif 520b3506946SBarry Smith 5216356e834SBarry Smith #undef __FUNCT__ 52253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private" 5238c34d3f5SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptions *PetscOptionsObject) 52453acd3b1SBarry Smith { 52553acd3b1SBarry Smith PetscErrorCode ierr; 5268c34d3f5SBarry Smith PetscOption last; 5276356e834SBarry Smith char option[256],value[1024],tmp[32]; 528330cf3c9SBarry Smith size_t j; 52953acd3b1SBarry Smith 53053acd3b1SBarry Smith PetscFunctionBegin; 53183355fc5SBarry Smith if (PetscOptionsObject->next) { 53283355fc5SBarry Smith if (!PetscOptionsObject->count) { 533a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS) 534f7b25cf6SBarry Smith ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr); 535b3506946SBarry Smith #else 536e55864a3SBarry Smith ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr); 537b3506946SBarry Smith #endif 538aee2cecaSBarry Smith } 539aee2cecaSBarry Smith } 5406356e834SBarry Smith 541e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr); 5426356e834SBarry Smith 543e26ddf31SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 544e55864a3SBarry Smith if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2; 5457a72a596SBarry Smith /* reset alreadyprinted flag */ 546e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = PETSC_FALSE; 547e55864a3SBarry Smith if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE; 548e55864a3SBarry Smith PetscOptionsObject->object = NULL; 54953acd3b1SBarry Smith 550e55864a3SBarry Smith while (PetscOptionsObject->next) { 551e55864a3SBarry Smith if (PetscOptionsObject->next->set) { 552e55864a3SBarry Smith if (PetscOptionsObject->prefix) { 55353acd3b1SBarry Smith ierr = PetscStrcpy(option,"-");CHKERRQ(ierr); 554e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr); 555e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr); 5566356e834SBarry Smith } else { 557e55864a3SBarry Smith ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr); 5586356e834SBarry Smith } 5596356e834SBarry Smith 560e55864a3SBarry Smith switch (PetscOptionsObject->next->type) { 5616356e834SBarry Smith case OPTION_HEAD: 5626356e834SBarry Smith break; 5636356e834SBarry Smith case OPTION_INT_ARRAY: 564e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]); 565e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 566e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]); 5676356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5686356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5696356e834SBarry Smith } 5706356e834SBarry Smith break; 5716356e834SBarry Smith case OPTION_INT: 572e55864a3SBarry Smith sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data); 5736356e834SBarry Smith break; 5746356e834SBarry Smith case OPTION_REAL: 575e55864a3SBarry Smith sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data); 5766356e834SBarry Smith break; 5776356e834SBarry Smith case OPTION_REAL_ARRAY: 578e55864a3SBarry Smith sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]); 579e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 580e55864a3SBarry Smith sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]); 5816356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5826356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5836356e834SBarry Smith } 5846356e834SBarry Smith break; 585050cccc3SHong Zhang case OPTION_SCALAR_ARRAY: 58695f3a755SJose E. Roman sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0])); 587050cccc3SHong Zhang for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 58895f3a755SJose E. Roman sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j])); 589050cccc3SHong Zhang ierr = PetscStrcat(value,",");CHKERRQ(ierr); 590050cccc3SHong Zhang ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 591050cccc3SHong Zhang } 592050cccc3SHong Zhang break; 5937781c08eSBarry Smith case OPTION_BOOL: 594e55864a3SBarry Smith sprintf(value,"%d",*(int*)PetscOptionsObject->next->data); 5956356e834SBarry Smith break; 5967781c08eSBarry Smith case OPTION_BOOL_ARRAY: 597e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]); 598e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 599e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]); 6001ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 6011ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 6021ae3d29cSBarry Smith } 6031ae3d29cSBarry Smith break; 604a264d7a6SBarry Smith case OPTION_FLIST: 6056991f827SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 6066991f827SBarry Smith break; 6071ae3d29cSBarry Smith case OPTION_ELIST: 608e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 6096356e834SBarry Smith break; 6101ae3d29cSBarry Smith case OPTION_STRING: 611e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 612d51da6bfSBarry Smith break; 6131ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 614e55864a3SBarry Smith sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]); 615e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 616e55864a3SBarry Smith sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]); 6171ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 6181ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 6191ae3d29cSBarry Smith } 6206356e834SBarry Smith break; 6216356e834SBarry Smith } 6226356e834SBarry Smith ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr); 6236356e834SBarry Smith } 6246991f827SBarry Smith if (PetscOptionsObject->next->type == OPTION_ELIST) { 6256991f827SBarry Smith ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr); 6266991f827SBarry Smith } 627e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr); 628e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr); 629e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr); 630e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr); 631c979a496SBarry Smith 63283355fc5SBarry Smith if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){ 63383355fc5SBarry Smith free(PetscOptionsObject->next->data); 634c979a496SBarry Smith } else { 63583355fc5SBarry Smith ierr = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr); 636c979a496SBarry Smith } 6377781c08eSBarry Smith 63883355fc5SBarry Smith last = PetscOptionsObject->next; 63983355fc5SBarry Smith PetscOptionsObject->next = PetscOptionsObject->next->next; 6406356e834SBarry Smith ierr = PetscFree(last);CHKERRQ(ierr); 6416356e834SBarry Smith } 642f59f755dSBarry Smith ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr); 643e55864a3SBarry Smith PetscOptionsObject->next = 0; 64453acd3b1SBarry Smith PetscFunctionReturn(0); 64553acd3b1SBarry Smith } 64653acd3b1SBarry Smith 64753acd3b1SBarry Smith #undef __FUNCT__ 648e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private" 64953acd3b1SBarry Smith /*@C 65053acd3b1SBarry Smith PetscOptionsEnum - Gets the enum value for a particular option in the database. 65153acd3b1SBarry Smith 6523f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 65353acd3b1SBarry Smith 65453acd3b1SBarry Smith Input Parameters: 65553acd3b1SBarry Smith + opt - option name 65653acd3b1SBarry Smith . text - short string that describes the option 65753acd3b1SBarry Smith . man - manual page with additional information on option 65853acd3b1SBarry Smith . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 6590fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6600fdccdaeSBarry Smith $ PetscOptionsEnum(..., obj->value,&object->value,...) or 6610fdccdaeSBarry Smith $ value = defaultvalue 6620fdccdaeSBarry Smith $ PetscOptionsEnum(..., value,&value,&flg); 6630fdccdaeSBarry Smith $ if (flg) { 66453acd3b1SBarry Smith 66553acd3b1SBarry Smith Output Parameter: 66653acd3b1SBarry Smith + value - the value to return 667b32e0204SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 66853acd3b1SBarry Smith 66953acd3b1SBarry Smith Level: beginner 67053acd3b1SBarry Smith 67153acd3b1SBarry Smith Concepts: options database 67253acd3b1SBarry Smith 67353acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 67453acd3b1SBarry Smith 67553acd3b1SBarry Smith list is usually something like PCASMTypes or some other predefined list of enum names 67653acd3b1SBarry Smith 67753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 678acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 679acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 68053acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 68153acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 682acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 683a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 68453acd3b1SBarry Smith @*/ 6858c34d3f5SBarry 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) 68653acd3b1SBarry Smith { 68753acd3b1SBarry Smith PetscErrorCode ierr; 68853acd3b1SBarry Smith PetscInt ntext = 0; 689aa5bb8c0SSatish Balay PetscInt tval; 690ace3abfcSBarry Smith PetscBool tflg; 69153acd3b1SBarry Smith 69253acd3b1SBarry Smith PetscFunctionBegin; 69353acd3b1SBarry Smith while (list[ntext++]) { 694e32f2f54SBarry Smith if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 69553acd3b1SBarry Smith } 696e32f2f54SBarry Smith if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 69753acd3b1SBarry Smith ntext -= 3; 698e55864a3SBarry Smith ierr = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr); 699aa5bb8c0SSatish Balay /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 700aa5bb8c0SSatish Balay if (tflg) *value = (PetscEnum)tval; 701aa5bb8c0SSatish Balay if (set) *set = tflg; 70253acd3b1SBarry Smith PetscFunctionReturn(0); 70353acd3b1SBarry Smith } 70453acd3b1SBarry Smith 705d3e47460SLisandro Dalcin #undef __FUNCT__ 706d3e47460SLisandro Dalcin #define __FUNCT__ "PetscOptionsEnumArray_Private" 707d3e47460SLisandro Dalcin /*@C 708d3e47460SLisandro Dalcin PetscOptionsEnumArray - Gets an array of enum values for a particular 709d3e47460SLisandro Dalcin option in the database. 710d3e47460SLisandro Dalcin 711d3e47460SLisandro Dalcin Logically Collective on the communicator passed in PetscOptionsBegin() 712d3e47460SLisandro Dalcin 713d3e47460SLisandro Dalcin Input Parameters: 714d3e47460SLisandro Dalcin + opt - the option one is seeking 715d3e47460SLisandro Dalcin . text - short string describing option 716d3e47460SLisandro Dalcin . man - manual page for option 717d3e47460SLisandro Dalcin - n - maximum number of values 718d3e47460SLisandro Dalcin 719d3e47460SLisandro Dalcin Output Parameter: 720d3e47460SLisandro Dalcin + value - location to copy values 721d3e47460SLisandro Dalcin . n - actual number of values found 722d3e47460SLisandro Dalcin - set - PETSC_TRUE if found, else PETSC_FALSE 723d3e47460SLisandro Dalcin 724d3e47460SLisandro Dalcin Level: beginner 725d3e47460SLisandro Dalcin 726d3e47460SLisandro Dalcin Notes: 727d3e47460SLisandro Dalcin The array must be passed as a comma separated list. 728d3e47460SLisandro Dalcin 729d3e47460SLisandro Dalcin There must be no intervening spaces between the values. 730d3e47460SLisandro Dalcin 731d3e47460SLisandro Dalcin Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 732d3e47460SLisandro Dalcin 733d3e47460SLisandro Dalcin Concepts: options database^array of enums 734d3e47460SLisandro Dalcin 735d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 736d3e47460SLisandro Dalcin PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 737d3e47460SLisandro Dalcin PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 738d3e47460SLisandro Dalcin PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 739d3e47460SLisandro Dalcin PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 740d3e47460SLisandro Dalcin PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray() 741d3e47460SLisandro Dalcin @*/ 742d3e47460SLisandro Dalcin PetscErrorCode PetscOptionsEnumArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool *set) 743d3e47460SLisandro Dalcin { 744d3e47460SLisandro Dalcin PetscInt i,nlist = 0; 745d3e47460SLisandro Dalcin PetscOption amsopt; 746d3e47460SLisandro Dalcin PetscErrorCode ierr; 747d3e47460SLisandro Dalcin 748d3e47460SLisandro Dalcin PetscFunctionBegin; 749d3e47460SLisandro Dalcin while (list[nlist++]) if (nlist > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 750d3e47460SLisandro Dalcin if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 751d3e47460SLisandro Dalcin nlist -= 3; /* drop enum name, prefix, and null termination */ 752d3e47460SLisandro Dalcin if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */ 753d3e47460SLisandro Dalcin PetscEnum *vals; 754d3e47460SLisandro Dalcin ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr); 755d3e47460SLisandro Dalcin ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr); 756d3e47460SLisandro Dalcin amsopt->nlist = nlist; 757d3e47460SLisandro Dalcin ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr); 758d3e47460SLisandro Dalcin amsopt->arraylength = *n; 759d3e47460SLisandro Dalcin vals = (PetscEnum*)amsopt->data; 760d3e47460SLisandro Dalcin for (i=0; i<*n; i++) vals[i] = value[i]; 761d3e47460SLisandro Dalcin } 762d3e47460SLisandro Dalcin ierr = PetscOptionsGetEnumArray(PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr); 763d3e47460SLisandro Dalcin if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 764d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr); 765d3e47460SLisandro Dalcin for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);} 766d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr); 767d3e47460SLisandro Dalcin for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);} 768d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr); 769d3e47460SLisandro Dalcin } 770d3e47460SLisandro Dalcin PetscFunctionReturn(0); 771d3e47460SLisandro Dalcin } 772d3e47460SLisandro Dalcin 77353acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/ 77453acd3b1SBarry Smith #undef __FUNCT__ 775e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private" 77653acd3b1SBarry Smith /*@C 77753acd3b1SBarry Smith PetscOptionsInt - Gets the integer value for a particular option in the database. 77853acd3b1SBarry Smith 7793f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 78053acd3b1SBarry Smith 78153acd3b1SBarry Smith Input Parameters: 78253acd3b1SBarry Smith + opt - option name 78353acd3b1SBarry Smith . text - short string that describes the option 78453acd3b1SBarry Smith . man - manual page with additional information on option 7850fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 7860fdccdaeSBarry Smith $ PetscOptionsInt(..., obj->value,&object->value,...) or 7870fdccdaeSBarry Smith $ value = defaultvalue 7880fdccdaeSBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 7890fdccdaeSBarry Smith $ if (flg) { 79053acd3b1SBarry Smith 79153acd3b1SBarry Smith Output Parameter: 79253acd3b1SBarry Smith + value - the integer value to return 79353acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 79453acd3b1SBarry Smith 79553acd3b1SBarry Smith Level: beginner 79653acd3b1SBarry Smith 79753acd3b1SBarry Smith Concepts: options database^has int 79853acd3b1SBarry Smith 79953acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 80053acd3b1SBarry Smith 80153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 802acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 803acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 80453acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 80553acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 806acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 807a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 80853acd3b1SBarry Smith @*/ 8098c34d3f5SBarry Smith PetscErrorCode PetscOptionsInt_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *set) 81053acd3b1SBarry Smith { 81153acd3b1SBarry Smith PetscErrorCode ierr; 8128c34d3f5SBarry Smith PetscOption amsopt; 81312655325SBarry Smith PetscBool wasset; 81453acd3b1SBarry Smith 81553acd3b1SBarry Smith PetscFunctionBegin; 816e55864a3SBarry Smith if (!PetscOptionsObject->count) { 817e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr); 8186356e834SBarry Smith ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr); 81912655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 8203e211508SBarry Smith 82112655325SBarry Smith ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,¤tvalue,&wasset);CHKERRQ(ierr); 8223e211508SBarry Smith if (wasset) { 82312655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 8243e211508SBarry Smith } 825af6d86caSBarry Smith } 826e55864a3SBarry Smith ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 827e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 8281a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 82953acd3b1SBarry Smith } 83053acd3b1SBarry Smith PetscFunctionReturn(0); 83153acd3b1SBarry Smith } 83253acd3b1SBarry Smith 83353acd3b1SBarry Smith #undef __FUNCT__ 834e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private" 83553acd3b1SBarry Smith /*@C 83653acd3b1SBarry Smith PetscOptionsString - Gets the string value for a particular option in the database. 83753acd3b1SBarry Smith 8383f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 83953acd3b1SBarry Smith 84053acd3b1SBarry Smith Input Parameters: 84153acd3b1SBarry Smith + opt - option name 84253acd3b1SBarry Smith . text - short string that describes the option 84353acd3b1SBarry Smith . man - manual page with additional information on option 8440fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 845bcbf2dc5SJed Brown - len - length of the result string including null terminator 84653acd3b1SBarry Smith 84753acd3b1SBarry Smith Output Parameter: 84853acd3b1SBarry Smith + value - the value to return 84953acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 85053acd3b1SBarry Smith 85153acd3b1SBarry Smith Level: beginner 85253acd3b1SBarry Smith 85353acd3b1SBarry Smith Concepts: options database^has int 85453acd3b1SBarry Smith 85553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 85653acd3b1SBarry Smith 8577fccdfe4SBarry 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). 8587fccdfe4SBarry Smith 85953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 860acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 861acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(), 86253acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 86353acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 864acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 865a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 86653acd3b1SBarry Smith @*/ 8678c34d3f5SBarry 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) 86853acd3b1SBarry Smith { 86953acd3b1SBarry Smith PetscErrorCode ierr; 8708c34d3f5SBarry Smith PetscOption amsopt; 87153acd3b1SBarry Smith 87253acd3b1SBarry Smith PetscFunctionBegin; 8731a1499c8SBarry Smith if (!PetscOptionsObject->count) { 8741a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 87564facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 8760fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 877af6d86caSBarry Smith } 878e55864a3SBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 879e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 8801a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 88153acd3b1SBarry Smith } 88253acd3b1SBarry Smith PetscFunctionReturn(0); 88353acd3b1SBarry Smith } 88453acd3b1SBarry Smith 88553acd3b1SBarry Smith #undef __FUNCT__ 886e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private" 88753acd3b1SBarry Smith /*@C 88853acd3b1SBarry Smith PetscOptionsReal - Gets the PetscReal value for a particular option in the database. 88953acd3b1SBarry Smith 8903f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 89153acd3b1SBarry Smith 89253acd3b1SBarry Smith Input Parameters: 89353acd3b1SBarry Smith + opt - option name 89453acd3b1SBarry Smith . text - short string that describes the option 89553acd3b1SBarry Smith . man - manual page with additional information on option 8960fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 8970fdccdaeSBarry Smith $ PetscOptionsReal(..., obj->value,&object->value,...) or 8980fdccdaeSBarry Smith $ value = defaultvalue 8990fdccdaeSBarry Smith $ PetscOptionsReal(..., value,&value,&flg); 9000fdccdaeSBarry Smith $ if (flg) { 90153acd3b1SBarry Smith 90253acd3b1SBarry Smith Output Parameter: 90353acd3b1SBarry Smith + value - the value to return 90453acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 90553acd3b1SBarry Smith 90653acd3b1SBarry Smith Level: beginner 90753acd3b1SBarry Smith 90853acd3b1SBarry Smith Concepts: options database^has int 90953acd3b1SBarry Smith 91053acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 91153acd3b1SBarry Smith 91253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 913acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 914acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 91553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 91653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 917acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 918a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 91953acd3b1SBarry Smith @*/ 9208c34d3f5SBarry Smith PetscErrorCode PetscOptionsReal_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 92153acd3b1SBarry Smith { 92253acd3b1SBarry Smith PetscErrorCode ierr; 9238c34d3f5SBarry Smith PetscOption amsopt; 92453acd3b1SBarry Smith 92553acd3b1SBarry Smith PetscFunctionBegin; 926e55864a3SBarry Smith if (!PetscOptionsObject->count) { 927e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr); 928538aa990SBarry Smith ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 929a297a907SKarl Rupp 9300fdccdaeSBarry Smith *(PetscReal*)amsopt->data = currentvalue; 931538aa990SBarry Smith } 9321a1499c8SBarry Smith ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 9331a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9341a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr); 93553acd3b1SBarry Smith } 93653acd3b1SBarry Smith PetscFunctionReturn(0); 93753acd3b1SBarry Smith } 93853acd3b1SBarry Smith 93953acd3b1SBarry Smith #undef __FUNCT__ 940e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private" 94153acd3b1SBarry Smith /*@C 94253acd3b1SBarry Smith PetscOptionsScalar - Gets the scalar value for a particular option in the database. 94353acd3b1SBarry Smith 9443f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 94553acd3b1SBarry Smith 94653acd3b1SBarry Smith Input Parameters: 94753acd3b1SBarry Smith + opt - option name 94853acd3b1SBarry Smith . text - short string that describes the option 94953acd3b1SBarry Smith . man - manual page with additional information on option 9500fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 9510fdccdaeSBarry Smith $ PetscOptionsScalar(..., obj->value,&object->value,...) or 9520fdccdaeSBarry Smith $ value = defaultvalue 9530fdccdaeSBarry Smith $ PetscOptionsScalar(..., value,&value,&flg); 9540fdccdaeSBarry Smith $ if (flg) { 9550fdccdaeSBarry Smith 95653acd3b1SBarry Smith 95753acd3b1SBarry Smith Output Parameter: 95853acd3b1SBarry Smith + value - the value to return 95953acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 96053acd3b1SBarry Smith 96153acd3b1SBarry Smith Level: beginner 96253acd3b1SBarry Smith 96353acd3b1SBarry Smith Concepts: options database^has int 96453acd3b1SBarry Smith 96553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 96653acd3b1SBarry Smith 96753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 968acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 969acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 97053acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 97153acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 972acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 973a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 97453acd3b1SBarry Smith @*/ 9758c34d3f5SBarry Smith PetscErrorCode PetscOptionsScalar_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 97653acd3b1SBarry Smith { 97753acd3b1SBarry Smith PetscErrorCode ierr; 97853acd3b1SBarry Smith 97953acd3b1SBarry Smith PetscFunctionBegin; 98053acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 9810fdccdaeSBarry Smith ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr); 98253acd3b1SBarry Smith #else 983e55864a3SBarry Smith ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 98453acd3b1SBarry Smith #endif 98553acd3b1SBarry Smith PetscFunctionReturn(0); 98653acd3b1SBarry Smith } 98753acd3b1SBarry Smith 98853acd3b1SBarry Smith #undef __FUNCT__ 989e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private" 99053acd3b1SBarry Smith /*@C 99190d69ab7SBarry 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 99290d69ab7SBarry Smith its value is set to false. 99353acd3b1SBarry Smith 9943f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 99553acd3b1SBarry Smith 99653acd3b1SBarry Smith Input Parameters: 99753acd3b1SBarry Smith + opt - option name 99853acd3b1SBarry Smith . text - short string that describes the option 99953acd3b1SBarry Smith - man - manual page with additional information on option 100053acd3b1SBarry Smith 100153acd3b1SBarry Smith Output Parameter: 100253acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 100353acd3b1SBarry Smith 100453acd3b1SBarry Smith Level: beginner 100553acd3b1SBarry Smith 100653acd3b1SBarry Smith Concepts: options database^has int 100753acd3b1SBarry Smith 100853acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 100953acd3b1SBarry Smith 101053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 1011acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1012acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 101353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 101453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1015acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1016a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 101753acd3b1SBarry Smith @*/ 10188c34d3f5SBarry Smith PetscErrorCode PetscOptionsName_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 101953acd3b1SBarry Smith { 102053acd3b1SBarry Smith PetscErrorCode ierr; 10218c34d3f5SBarry Smith PetscOption amsopt; 102253acd3b1SBarry Smith 102353acd3b1SBarry Smith PetscFunctionBegin; 1024e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1025e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1026ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1027a297a907SKarl Rupp 1028ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 10291ae3d29cSBarry Smith } 1030e55864a3SBarry Smith ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr); 1031e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1032e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 103353acd3b1SBarry Smith } 103453acd3b1SBarry Smith PetscFunctionReturn(0); 103553acd3b1SBarry Smith } 103653acd3b1SBarry Smith 103753acd3b1SBarry Smith #undef __FUNCT__ 1038e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private" 103953acd3b1SBarry Smith /*@C 1040a264d7a6SBarry Smith PetscOptionsFList - Puts a list of option values that a single one may be selected from 104153acd3b1SBarry Smith 10423f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 104353acd3b1SBarry Smith 104453acd3b1SBarry Smith Input Parameters: 104553acd3b1SBarry Smith + opt - option name 104653acd3b1SBarry Smith . text - short string that describes the option 104753acd3b1SBarry Smith . man - manual page with additional information on option 104853acd3b1SBarry Smith . list - the possible choices 10490fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10500fdccdaeSBarry Smith $ PetscOptionsFlist(..., obj->value,value,len,&flg); 10510fdccdaeSBarry Smith $ if (flg) { 10523cc1e11dSBarry Smith - len - the length of the character array value 105353acd3b1SBarry Smith 105453acd3b1SBarry Smith Output Parameter: 105553acd3b1SBarry Smith + value - the value to return 105653acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 105753acd3b1SBarry Smith 105853acd3b1SBarry Smith Level: intermediate 105953acd3b1SBarry Smith 106053acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 106153acd3b1SBarry Smith 106253acd3b1SBarry Smith See PetscOptionsEList() for when the choices are given in a string array 106353acd3b1SBarry Smith 106453acd3b1SBarry Smith To get a listing of all currently specified options, 106588c29154SBarry Smith see PetscOptionsView() or PetscOptionsGetAll() 106653acd3b1SBarry Smith 1067eabe10d7SBarry Smith Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list 1068eabe10d7SBarry Smith 106953acd3b1SBarry Smith Concepts: options database^list 107053acd3b1SBarry Smith 107153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1072acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 107353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 107453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1075acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1076a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum() 107753acd3b1SBarry Smith @*/ 10788c34d3f5SBarry 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) 107953acd3b1SBarry Smith { 108053acd3b1SBarry Smith PetscErrorCode ierr; 10818c34d3f5SBarry Smith PetscOption amsopt; 108253acd3b1SBarry Smith 108353acd3b1SBarry Smith PetscFunctionBegin; 10841a1499c8SBarry Smith if (!PetscOptionsObject->count) { 10851a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr); 108664facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 10870fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 10883cc1e11dSBarry Smith amsopt->flist = list; 10893cc1e11dSBarry Smith } 10901a1499c8SBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 10911a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10921a1499c8SBarry Smith ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr); 109353acd3b1SBarry Smith } 109453acd3b1SBarry Smith PetscFunctionReturn(0); 109553acd3b1SBarry Smith } 109653acd3b1SBarry Smith 109753acd3b1SBarry Smith #undef __FUNCT__ 1098e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private" 109953acd3b1SBarry Smith /*@C 110053acd3b1SBarry Smith PetscOptionsEList - Puts a list of option values that a single one may be selected from 110153acd3b1SBarry Smith 11023f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 110353acd3b1SBarry Smith 110453acd3b1SBarry Smith Input Parameters: 110553acd3b1SBarry Smith + opt - option name 110653acd3b1SBarry Smith . ltext - short string that describes the option 110753acd3b1SBarry Smith . man - manual page with additional information on option 1108a264d7a6SBarry Smith . list - the possible choices (one of these must be selected, anything else is invalid) 110953acd3b1SBarry Smith . ntext - number of choices 11100fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 11110fdccdaeSBarry Smith $ PetscOptionsElist(..., obj->value,&value,&flg); 11120fdccdaeSBarry Smith $ if (flg) { 11130fdccdaeSBarry Smith 111453acd3b1SBarry Smith 111553acd3b1SBarry Smith Output Parameter: 111653acd3b1SBarry Smith + value - the index of the value to return 111753acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 111853acd3b1SBarry Smith 111953acd3b1SBarry Smith Level: intermediate 112053acd3b1SBarry Smith 112153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 112253acd3b1SBarry Smith 1123a264d7a6SBarry Smith See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 112453acd3b1SBarry Smith 112553acd3b1SBarry Smith Concepts: options database^list 112653acd3b1SBarry Smith 112753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1128acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 112953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 113053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1131acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1132a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEnum() 113353acd3b1SBarry Smith @*/ 11348c34d3f5SBarry 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) 113553acd3b1SBarry Smith { 113653acd3b1SBarry Smith PetscErrorCode ierr; 113753acd3b1SBarry Smith PetscInt i; 11388c34d3f5SBarry Smith PetscOption amsopt; 113953acd3b1SBarry Smith 114053acd3b1SBarry Smith PetscFunctionBegin; 11411a1499c8SBarry Smith if (!PetscOptionsObject->count) { 11421a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr); 114364facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 11440fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 11456991f827SBarry Smith ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr); 11461ae3d29cSBarry Smith amsopt->nlist = ntext; 11471ae3d29cSBarry Smith } 11481a1499c8SBarry Smith ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr); 11491a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 11501a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr); 115153acd3b1SBarry Smith for (i=0; i<ntext; i++) { 1152e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr); 115353acd3b1SBarry Smith } 1154e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr); 115553acd3b1SBarry Smith } 115653acd3b1SBarry Smith PetscFunctionReturn(0); 115753acd3b1SBarry Smith } 115853acd3b1SBarry Smith 115953acd3b1SBarry Smith #undef __FUNCT__ 1160e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private" 116153acd3b1SBarry Smith /*@C 1162acfcf0e5SJed Brown PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 1163d5649816SBarry Smith which at most a single value can be true. 116453acd3b1SBarry Smith 11653f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 116653acd3b1SBarry Smith 116753acd3b1SBarry Smith Input Parameters: 116853acd3b1SBarry Smith + opt - option name 116953acd3b1SBarry Smith . text - short string that describes the option 117053acd3b1SBarry Smith - man - manual page with additional information on option 117153acd3b1SBarry Smith 117253acd3b1SBarry Smith Output Parameter: 117353acd3b1SBarry Smith . flg - whether that option was set or not 117453acd3b1SBarry Smith 117553acd3b1SBarry Smith Level: intermediate 117653acd3b1SBarry Smith 117753acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 117853acd3b1SBarry Smith 1179acfcf0e5SJed Brown Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd() 118053acd3b1SBarry Smith 118153acd3b1SBarry Smith Concepts: options database^logical group 118253acd3b1SBarry Smith 118353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1184acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 118553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 118653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1187acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1188a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 118953acd3b1SBarry Smith @*/ 11908c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 119153acd3b1SBarry Smith { 119253acd3b1SBarry Smith PetscErrorCode ierr; 11938c34d3f5SBarry Smith PetscOption amsopt; 119453acd3b1SBarry Smith 119553acd3b1SBarry Smith PetscFunctionBegin; 1196e55864a3SBarry Smith if (!PetscOptionsObject->count) { 119783355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1198ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1199a297a907SKarl Rupp 1200ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 12011ae3d29cSBarry Smith } 120268b16fdaSBarry Smith *flg = PETSC_FALSE; 1203e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1204e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1205e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," Pick at most one of -------------\n");CHKERRQ(ierr); 1206e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 120753acd3b1SBarry Smith } 120853acd3b1SBarry Smith PetscFunctionReturn(0); 120953acd3b1SBarry Smith } 121053acd3b1SBarry Smith 121153acd3b1SBarry Smith #undef __FUNCT__ 1212e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private" 121353acd3b1SBarry Smith /*@C 1214acfcf0e5SJed Brown PetscOptionsBoolGroup - One in a series of logical queries on the options database for 1215d5649816SBarry Smith which at most a single value can be true. 121653acd3b1SBarry Smith 12173f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 121853acd3b1SBarry Smith 121953acd3b1SBarry Smith Input Parameters: 122053acd3b1SBarry Smith + opt - option name 122153acd3b1SBarry Smith . text - short string that describes the option 122253acd3b1SBarry Smith - man - manual page with additional information on option 122353acd3b1SBarry Smith 122453acd3b1SBarry Smith Output Parameter: 122553acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 122653acd3b1SBarry Smith 122753acd3b1SBarry Smith Level: intermediate 122853acd3b1SBarry Smith 122953acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 123053acd3b1SBarry Smith 1231acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd() 123253acd3b1SBarry Smith 123353acd3b1SBarry Smith Concepts: options database^logical group 123453acd3b1SBarry Smith 123553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1236acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 123753acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 123853acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1239acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1240a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 124153acd3b1SBarry Smith @*/ 12428c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 124353acd3b1SBarry Smith { 124453acd3b1SBarry Smith PetscErrorCode ierr; 12458c34d3f5SBarry Smith PetscOption amsopt; 124653acd3b1SBarry Smith 124753acd3b1SBarry Smith PetscFunctionBegin; 1248e55864a3SBarry Smith if (!PetscOptionsObject->count) { 124983355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1250ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1251a297a907SKarl Rupp 1252ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 12531ae3d29cSBarry Smith } 125417326d04SJed Brown *flg = PETSC_FALSE; 1255e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1256e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1257e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 125853acd3b1SBarry Smith } 125953acd3b1SBarry Smith PetscFunctionReturn(0); 126053acd3b1SBarry Smith } 126153acd3b1SBarry Smith 126253acd3b1SBarry Smith #undef __FUNCT__ 1263e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private" 126453acd3b1SBarry Smith /*@C 1265acfcf0e5SJed Brown PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 1266d5649816SBarry Smith which at most a single value can be true. 126753acd3b1SBarry Smith 12683f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 126953acd3b1SBarry Smith 127053acd3b1SBarry Smith Input Parameters: 127153acd3b1SBarry Smith + opt - option name 127253acd3b1SBarry Smith . text - short string that describes the option 127353acd3b1SBarry Smith - man - manual page with additional information on option 127453acd3b1SBarry Smith 127553acd3b1SBarry Smith Output Parameter: 127653acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 127753acd3b1SBarry Smith 127853acd3b1SBarry Smith Level: intermediate 127953acd3b1SBarry Smith 128053acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 128153acd3b1SBarry Smith 1282acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() 128353acd3b1SBarry Smith 128453acd3b1SBarry Smith Concepts: options database^logical group 128553acd3b1SBarry Smith 128653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1287acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 128853acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 128953acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1290acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1291a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 129253acd3b1SBarry Smith @*/ 12938c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 129453acd3b1SBarry Smith { 129553acd3b1SBarry Smith PetscErrorCode ierr; 12968c34d3f5SBarry Smith PetscOption amsopt; 129753acd3b1SBarry Smith 129853acd3b1SBarry Smith PetscFunctionBegin; 1299e55864a3SBarry Smith if (!PetscOptionsObject->count) { 130083355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1301ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1302a297a907SKarl Rupp 1303ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 13041ae3d29cSBarry Smith } 130517326d04SJed Brown *flg = PETSC_FALSE; 1306e55864a3SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1307e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1308e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 130953acd3b1SBarry Smith } 131053acd3b1SBarry Smith PetscFunctionReturn(0); 131153acd3b1SBarry Smith } 131253acd3b1SBarry Smith 131353acd3b1SBarry Smith #undef __FUNCT__ 1314e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private" 131553acd3b1SBarry Smith /*@C 1316acfcf0e5SJed Brown PetscOptionsBool - Determines if a particular option is in the database with a true or false 131753acd3b1SBarry Smith 13183f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 131953acd3b1SBarry Smith 132053acd3b1SBarry Smith Input Parameters: 132153acd3b1SBarry Smith + opt - option name 132253acd3b1SBarry Smith . text - short string that describes the option 1323868c398cSBarry Smith . man - manual page with additional information on option 132494ae4db5SBarry Smith - currentvalue - the current value 132553acd3b1SBarry Smith 132653acd3b1SBarry Smith Output Parameter: 132753acd3b1SBarry Smith . flg - PETSC_TRUE or PETSC_FALSE 132853acd3b1SBarry Smith . set - PETSC_TRUE if found, else PETSC_FALSE 132953acd3b1SBarry Smith 133053acd3b1SBarry Smith Level: beginner 133153acd3b1SBarry Smith 133253acd3b1SBarry Smith Concepts: options database^logical 133353acd3b1SBarry Smith 133453acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 133553acd3b1SBarry Smith 133653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 1337acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1338acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 133953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 134053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1341acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1342a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 134353acd3b1SBarry Smith @*/ 13448c34d3f5SBarry Smith PetscErrorCode PetscOptionsBool_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 134553acd3b1SBarry Smith { 134653acd3b1SBarry Smith PetscErrorCode ierr; 1347ace3abfcSBarry Smith PetscBool iset; 13488c34d3f5SBarry Smith PetscOption amsopt; 134953acd3b1SBarry Smith 135053acd3b1SBarry Smith PetscFunctionBegin; 1351e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1352e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1353ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1354a297a907SKarl Rupp 135594ae4db5SBarry Smith *(PetscBool*)amsopt->data = currentvalue; 1356af6d86caSBarry Smith } 13571a1499c8SBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr); 135853acd3b1SBarry Smith if (set) *set = iset; 13591a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 136094ae4db5SBarry Smith const char *v = PetscBools[currentvalue]; 13611a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr); 136253acd3b1SBarry Smith } 136353acd3b1SBarry Smith PetscFunctionReturn(0); 136453acd3b1SBarry Smith } 136553acd3b1SBarry Smith 136653acd3b1SBarry Smith #undef __FUNCT__ 1367e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private" 136853acd3b1SBarry Smith /*@C 136953acd3b1SBarry Smith PetscOptionsRealArray - Gets an array of double values for a particular 137053acd3b1SBarry Smith option in the database. The values must be separated with commas with 137153acd3b1SBarry Smith no intervening spaces. 137253acd3b1SBarry Smith 13733f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 137453acd3b1SBarry Smith 137553acd3b1SBarry Smith Input Parameters: 137653acd3b1SBarry Smith + opt - the option one is seeking 137753acd3b1SBarry Smith . text - short string describing option 137853acd3b1SBarry Smith . man - manual page for option 137953acd3b1SBarry Smith - nmax - maximum number of values 138053acd3b1SBarry Smith 138153acd3b1SBarry Smith Output Parameter: 138253acd3b1SBarry Smith + value - location to copy values 138353acd3b1SBarry Smith . nmax - actual number of values found 138453acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 138553acd3b1SBarry Smith 138653acd3b1SBarry Smith Level: beginner 138753acd3b1SBarry Smith 138853acd3b1SBarry Smith Notes: 138953acd3b1SBarry Smith The user should pass in an array of doubles 139053acd3b1SBarry Smith 139153acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 139253acd3b1SBarry Smith 139353acd3b1SBarry Smith Concepts: options database^array of strings 139453acd3b1SBarry Smith 139553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1396acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 139753acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 139853acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1399acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1400a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 140153acd3b1SBarry Smith @*/ 14028c34d3f5SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 140353acd3b1SBarry Smith { 140453acd3b1SBarry Smith PetscErrorCode ierr; 140553acd3b1SBarry Smith PetscInt i; 14068c34d3f5SBarry Smith PetscOption amsopt; 140753acd3b1SBarry Smith 140853acd3b1SBarry Smith PetscFunctionBegin; 1409e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1410e26ddf31SBarry Smith PetscReal *vals; 1411e26ddf31SBarry Smith 1412e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr); 1413e55864a3SBarry Smith ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 1414e26ddf31SBarry Smith vals = (PetscReal*)amsopt->data; 1415e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1416e26ddf31SBarry Smith amsopt->arraylength = *n; 1417e26ddf31SBarry Smith } 1418e55864a3SBarry Smith ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1419e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1420a519f713SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr); 142153acd3b1SBarry Smith for (i=1; i<*n; i++) { 1422e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr); 142353acd3b1SBarry Smith } 1424e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 142553acd3b1SBarry Smith } 142653acd3b1SBarry Smith PetscFunctionReturn(0); 142753acd3b1SBarry Smith } 142853acd3b1SBarry Smith 1429050cccc3SHong Zhang #undef __FUNCT__ 1430050cccc3SHong Zhang #define __FUNCT__ "PetscOptionsScalarArray_Private" 1431050cccc3SHong Zhang /*@C 1432050cccc3SHong Zhang PetscOptionsScalarArray - Gets an array of Scalar values for a particular 1433050cccc3SHong Zhang option in the database. The values must be separated with commas with 1434050cccc3SHong Zhang no intervening spaces. 1435050cccc3SHong Zhang 1436050cccc3SHong Zhang Logically Collective on the communicator passed in PetscOptionsBegin() 1437050cccc3SHong Zhang 1438050cccc3SHong Zhang Input Parameters: 1439050cccc3SHong Zhang + opt - the option one is seeking 1440050cccc3SHong Zhang . text - short string describing option 1441050cccc3SHong Zhang . man - manual page for option 1442050cccc3SHong Zhang - nmax - maximum number of values 1443050cccc3SHong Zhang 1444050cccc3SHong Zhang Output Parameter: 1445050cccc3SHong Zhang + value - location to copy values 1446050cccc3SHong Zhang . nmax - actual number of values found 1447050cccc3SHong Zhang - set - PETSC_TRUE if found, else PETSC_FALSE 1448050cccc3SHong Zhang 1449050cccc3SHong Zhang Level: beginner 1450050cccc3SHong Zhang 1451050cccc3SHong Zhang Notes: 1452050cccc3SHong Zhang The user should pass in an array of doubles 1453050cccc3SHong Zhang 1454050cccc3SHong Zhang Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1455050cccc3SHong Zhang 1456050cccc3SHong Zhang Concepts: options database^array of strings 1457050cccc3SHong Zhang 1458050cccc3SHong Zhang .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1459050cccc3SHong Zhang PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1460050cccc3SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1461050cccc3SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1462050cccc3SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1463050cccc3SHong Zhang PetscOptionsFList(), PetscOptionsEList() 1464050cccc3SHong Zhang @*/ 1465050cccc3SHong Zhang PetscErrorCode PetscOptionsScalarArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool *set) 1466050cccc3SHong Zhang { 1467050cccc3SHong Zhang PetscErrorCode ierr; 1468050cccc3SHong Zhang PetscInt i; 1469050cccc3SHong Zhang PetscOption amsopt; 1470050cccc3SHong Zhang 1471050cccc3SHong Zhang PetscFunctionBegin; 1472050cccc3SHong Zhang if (!PetscOptionsObject->count) { 1473050cccc3SHong Zhang PetscScalar *vals; 1474050cccc3SHong Zhang 1475050cccc3SHong Zhang ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr); 1476050cccc3SHong Zhang ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr); 1477050cccc3SHong Zhang vals = (PetscScalar*)amsopt->data; 1478050cccc3SHong Zhang for (i=0; i<*n; i++) vals[i] = value[i]; 1479050cccc3SHong Zhang amsopt->arraylength = *n; 1480050cccc3SHong Zhang } 1481050cccc3SHong Zhang ierr = PetscOptionsGetScalarArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1482050cccc3SHong Zhang if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 148395f3a755SJose E. Roman ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g+%gi",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)PetscRealPart(value[0]),(double)PetscImaginaryPart(value[0]));CHKERRQ(ierr); 1484050cccc3SHong Zhang for (i=1; i<*n; i++) { 148595f3a755SJose E. Roman ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr); 1486050cccc3SHong Zhang } 1487050cccc3SHong Zhang ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 1488050cccc3SHong Zhang } 1489050cccc3SHong Zhang PetscFunctionReturn(0); 1490050cccc3SHong Zhang } 149153acd3b1SBarry Smith 149253acd3b1SBarry Smith #undef __FUNCT__ 1493e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private" 149453acd3b1SBarry Smith /*@C 149553acd3b1SBarry Smith PetscOptionsIntArray - Gets an array of integers for a particular 1496b32a342fSShri Abhyankar option in the database. 149753acd3b1SBarry Smith 14983f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 149953acd3b1SBarry Smith 150053acd3b1SBarry Smith Input Parameters: 150153acd3b1SBarry Smith + opt - the option one is seeking 150253acd3b1SBarry Smith . text - short string describing option 150353acd3b1SBarry Smith . man - manual page for option 1504f8a50e2bSBarry Smith - n - maximum number of values 150553acd3b1SBarry Smith 150653acd3b1SBarry Smith Output Parameter: 150753acd3b1SBarry Smith + value - location to copy values 1508f8a50e2bSBarry Smith . n - actual number of values found 150953acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 151053acd3b1SBarry Smith 151153acd3b1SBarry Smith Level: beginner 151253acd3b1SBarry Smith 151353acd3b1SBarry Smith Notes: 1514b32a342fSShri Abhyankar The array can be passed as 1515*bebe2cf6SSatish Balay a comma separated list: 0,1,2,3,4,5,6,7 15160fd488f5SShri Abhyankar a range (start-end+1): 0-8 15170fd488f5SShri Abhyankar a range with given increment (start-end+1:inc): 0-7:2 1518*bebe2cf6SSatish Balay a combination of values and ranges separated by commas: 0,1-8,8-15:2 1519b32a342fSShri Abhyankar 1520b32a342fSShri Abhyankar There must be no intervening spaces between the values. 152153acd3b1SBarry Smith 152253acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 152353acd3b1SBarry Smith 1524b32a342fSShri Abhyankar Concepts: options database^array of ints 152553acd3b1SBarry Smith 152653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1527acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 152853acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 152953acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1530acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1531a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray() 153253acd3b1SBarry Smith @*/ 15338c34d3f5SBarry Smith PetscErrorCode PetscOptionsIntArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 153453acd3b1SBarry Smith { 153553acd3b1SBarry Smith PetscErrorCode ierr; 153653acd3b1SBarry Smith PetscInt i; 15378c34d3f5SBarry Smith PetscOption amsopt; 153853acd3b1SBarry Smith 153953acd3b1SBarry Smith PetscFunctionBegin; 1540e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1541e26ddf31SBarry Smith PetscInt *vals; 1542e26ddf31SBarry Smith 1543e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr); 1544854ce69bSBarry Smith ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr); 1545e26ddf31SBarry Smith vals = (PetscInt*)amsopt->data; 1546e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1547e26ddf31SBarry Smith amsopt->arraylength = *n; 1548e26ddf31SBarry Smith } 1549e55864a3SBarry Smith ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1550e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1551e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 155253acd3b1SBarry Smith for (i=1; i<*n; i++) { 1553e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 155453acd3b1SBarry Smith } 1555e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 155653acd3b1SBarry Smith } 155753acd3b1SBarry Smith PetscFunctionReturn(0); 155853acd3b1SBarry Smith } 155953acd3b1SBarry Smith 156053acd3b1SBarry Smith #undef __FUNCT__ 1561e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private" 156253acd3b1SBarry Smith /*@C 156353acd3b1SBarry Smith PetscOptionsStringArray - Gets an array of string values for a particular 156453acd3b1SBarry Smith option in the database. The values must be separated with commas with 156553acd3b1SBarry Smith no intervening spaces. 156653acd3b1SBarry Smith 15673f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 156853acd3b1SBarry Smith 156953acd3b1SBarry Smith Input Parameters: 157053acd3b1SBarry Smith + opt - the option one is seeking 157153acd3b1SBarry Smith . text - short string describing option 157253acd3b1SBarry Smith . man - manual page for option 157353acd3b1SBarry Smith - nmax - maximum number of strings 157453acd3b1SBarry Smith 157553acd3b1SBarry Smith Output Parameter: 157653acd3b1SBarry Smith + value - location to copy strings 157753acd3b1SBarry Smith . nmax - actual number of strings found 157853acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 157953acd3b1SBarry Smith 158053acd3b1SBarry Smith Level: beginner 158153acd3b1SBarry Smith 158253acd3b1SBarry Smith Notes: 158353acd3b1SBarry Smith The user should pass in an array of pointers to char, to hold all the 158453acd3b1SBarry Smith strings returned by this function. 158553acd3b1SBarry Smith 158653acd3b1SBarry Smith The user is responsible for deallocating the strings that are 158753acd3b1SBarry Smith returned. The Fortran interface for this routine is not supported. 158853acd3b1SBarry Smith 158953acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 159053acd3b1SBarry Smith 159153acd3b1SBarry Smith Concepts: options database^array of strings 159253acd3b1SBarry Smith 159353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1594acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 159553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 159653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1597acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1598a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 159953acd3b1SBarry Smith @*/ 16008c34d3f5SBarry Smith PetscErrorCode PetscOptionsStringArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 160153acd3b1SBarry Smith { 160253acd3b1SBarry Smith PetscErrorCode ierr; 16038c34d3f5SBarry Smith PetscOption amsopt; 160453acd3b1SBarry Smith 160553acd3b1SBarry Smith PetscFunctionBegin; 1606e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1607e55864a3SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr); 1608854ce69bSBarry Smith ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr); 1609a297a907SKarl Rupp 16101ae3d29cSBarry Smith amsopt->arraylength = *nmax; 16111ae3d29cSBarry Smith } 1612e55864a3SBarry Smith ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr); 1613e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1614e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 161553acd3b1SBarry Smith } 161653acd3b1SBarry Smith PetscFunctionReturn(0); 161753acd3b1SBarry Smith } 161853acd3b1SBarry Smith 1619e2446a98SMatthew Knepley #undef __FUNCT__ 1620e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private" 1621e2446a98SMatthew Knepley /*@C 1622acfcf0e5SJed Brown PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 1623e2446a98SMatthew Knepley option in the database. The values must be separated with commas with 1624e2446a98SMatthew Knepley no intervening spaces. 1625e2446a98SMatthew Knepley 16263f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 1627e2446a98SMatthew Knepley 1628e2446a98SMatthew Knepley Input Parameters: 1629e2446a98SMatthew Knepley + opt - the option one is seeking 1630e2446a98SMatthew Knepley . text - short string describing option 1631e2446a98SMatthew Knepley . man - manual page for option 1632e2446a98SMatthew Knepley - nmax - maximum number of values 1633e2446a98SMatthew Knepley 1634e2446a98SMatthew Knepley Output Parameter: 1635e2446a98SMatthew Knepley + value - location to copy values 1636e2446a98SMatthew Knepley . nmax - actual number of values found 1637e2446a98SMatthew Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 1638e2446a98SMatthew Knepley 1639e2446a98SMatthew Knepley Level: beginner 1640e2446a98SMatthew Knepley 1641e2446a98SMatthew Knepley Notes: 1642e2446a98SMatthew Knepley The user should pass in an array of doubles 1643e2446a98SMatthew Knepley 1644e2446a98SMatthew Knepley Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1645e2446a98SMatthew Knepley 1646e2446a98SMatthew Knepley Concepts: options database^array of strings 1647e2446a98SMatthew Knepley 1648e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1649acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1650e2446a98SMatthew Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1651e2446a98SMatthew Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1652acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1653a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 1654e2446a98SMatthew Knepley @*/ 16558c34d3f5SBarry Smith PetscErrorCode PetscOptionsBoolArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 1656e2446a98SMatthew Knepley { 1657e2446a98SMatthew Knepley PetscErrorCode ierr; 1658e2446a98SMatthew Knepley PetscInt i; 16598c34d3f5SBarry Smith PetscOption amsopt; 1660e2446a98SMatthew Knepley 1661e2446a98SMatthew Knepley PetscFunctionBegin; 1662e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1663ace3abfcSBarry Smith PetscBool *vals; 16641ae3d29cSBarry Smith 166583355fc5SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr); 16661a1499c8SBarry Smith ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr); 1667ace3abfcSBarry Smith vals = (PetscBool*)amsopt->data; 16681ae3d29cSBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 16691ae3d29cSBarry Smith amsopt->arraylength = *n; 16701ae3d29cSBarry Smith } 1671e55864a3SBarry Smith ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1672e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1673e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 1674e2446a98SMatthew Knepley for (i=1; i<*n; i++) { 1675e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 1676e2446a98SMatthew Knepley } 1677e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 1678e2446a98SMatthew Knepley } 1679e2446a98SMatthew Knepley PetscFunctionReturn(0); 1680e2446a98SMatthew Knepley } 1681e2446a98SMatthew Knepley 16828cc676e6SMatthew G Knepley #undef __FUNCT__ 1683e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private" 16848cc676e6SMatthew G Knepley /*@C 1685d1da0b69SBarry Smith PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user 16868cc676e6SMatthew G Knepley 16878cc676e6SMatthew G Knepley Logically Collective on the communicator passed in PetscOptionsBegin() 16888cc676e6SMatthew G Knepley 16898cc676e6SMatthew G Knepley Input Parameters: 16908cc676e6SMatthew G Knepley + opt - option name 16918cc676e6SMatthew G Knepley . text - short string that describes the option 16928cc676e6SMatthew G Knepley - man - manual page with additional information on option 16938cc676e6SMatthew G Knepley 16948cc676e6SMatthew G Knepley Output Parameter: 16958cc676e6SMatthew G Knepley + viewer - the viewer 16968cc676e6SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 16978cc676e6SMatthew G Knepley 16988cc676e6SMatthew G Knepley Level: beginner 16998cc676e6SMatthew G Knepley 17008cc676e6SMatthew G Knepley Concepts: options database^has int 17018cc676e6SMatthew G Knepley 17028cc676e6SMatthew G Knepley Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 17038cc676e6SMatthew G Knepley 17045a7113b9SPatrick Sanan See PetscOptionsGetViewer() for the format of the supplied viewer and its options 17058cc676e6SMatthew G Knepley 17068cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 17078cc676e6SMatthew G Knepley PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 17088cc676e6SMatthew G Knepley PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 17098cc676e6SMatthew G Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 17108cc676e6SMatthew G Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 17118cc676e6SMatthew G Knepley PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1712a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 17138cc676e6SMatthew G Knepley @*/ 17148c34d3f5SBarry Smith PetscErrorCode PetscOptionsViewer_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 17158cc676e6SMatthew G Knepley { 17168cc676e6SMatthew G Knepley PetscErrorCode ierr; 17178c34d3f5SBarry Smith PetscOption amsopt; 17188cc676e6SMatthew G Knepley 17198cc676e6SMatthew G Knepley PetscFunctionBegin; 17201a1499c8SBarry Smith if (!PetscOptionsObject->count) { 17211a1499c8SBarry Smith ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 172264facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 17235b02f95dSBarry Smith ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr); 17248cc676e6SMatthew G Knepley } 1725e55864a3SBarry Smith ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr); 1726e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1727e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr); 17288cc676e6SMatthew G Knepley } 17298cc676e6SMatthew G Knepley PetscFunctionReturn(0); 17308cc676e6SMatthew G Knepley } 17318cc676e6SMatthew G Knepley 173253acd3b1SBarry Smith 173353acd3b1SBarry Smith #undef __FUNCT__ 173453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead" 173553acd3b1SBarry Smith /*@C 1736b52f573bSBarry Smith PetscOptionsHead - Puts a heading before listing any more published options. Used, for example, 173753acd3b1SBarry Smith in KSPSetFromOptions_GMRES(). 173853acd3b1SBarry Smith 17393f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 174053acd3b1SBarry Smith 174153acd3b1SBarry Smith Input Parameter: 174253acd3b1SBarry Smith . head - the heading text 174353acd3b1SBarry Smith 174453acd3b1SBarry Smith 174553acd3b1SBarry Smith Level: intermediate 174653acd3b1SBarry Smith 174753acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 174853acd3b1SBarry Smith 1749b52f573bSBarry Smith Can be followed by a call to PetscOptionsTail() in the same function. 175053acd3b1SBarry Smith 175153acd3b1SBarry Smith Concepts: options database^subheading 175253acd3b1SBarry Smith 175353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 1754acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 175553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 175653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1757acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1758a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 175953acd3b1SBarry Smith @*/ 17608c34d3f5SBarry Smith PetscErrorCode PetscOptionsHead(PetscOptions *PetscOptionsObject,const char head[]) 176153acd3b1SBarry Smith { 176253acd3b1SBarry Smith PetscErrorCode ierr; 176353acd3b1SBarry Smith 176453acd3b1SBarry Smith PetscFunctionBegin; 1765e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1766e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s\n",head);CHKERRQ(ierr); 176753acd3b1SBarry Smith } 176853acd3b1SBarry Smith PetscFunctionReturn(0); 176953acd3b1SBarry Smith } 177053acd3b1SBarry Smith 177153acd3b1SBarry Smith 177253acd3b1SBarry Smith 177353acd3b1SBarry Smith 177453acd3b1SBarry Smith 177553acd3b1SBarry Smith 1776