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 /* 2353acd3b1SBarry Smith Handles setting up the data structure in a call to PetscOptionsBegin() 2453acd3b1SBarry Smith */ 254416b707SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[]) 2653acd3b1SBarry Smith { 2753acd3b1SBarry Smith PetscErrorCode ierr; 2853acd3b1SBarry Smith 2953acd3b1SBarry Smith PetscFunctionBegin; 300eb63584SBarry Smith if (!PetscOptionsObject->alreadyprinted) { 319de0f6ecSBarry Smith if (!PetscOptionsHelpPrintedSingleton) { 329de0f6ecSBarry Smith ierr = PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);CHKERRQ(ierr); 339de0f6ecSBarry Smith } 349de0f6ecSBarry Smith ierr = PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,prefix,title,&PetscOptionsObject->alreadyprinted);CHKERRQ(ierr); 350eb63584SBarry Smith } 36e55864a3SBarry Smith PetscOptionsObject->next = 0; 37e55864a3SBarry Smith PetscOptionsObject->comm = comm; 38e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_FALSE; 39a297a907SKarl Rupp 40e55864a3SBarry Smith ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr); 41e55864a3SBarry Smith ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr); 4253acd3b1SBarry Smith 43c5929fdfSBarry Smith ierr = PetscOptionsHasName(PetscOptionsObject->options,NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr); 44e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) { 45e55864a3SBarry Smith if (!PetscOptionsObject->alreadyprinted) { 4653acd3b1SBarry Smith ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr); 4753acd3b1SBarry Smith } 4861b37b28SSatish Balay } 4953acd3b1SBarry Smith PetscFunctionReturn(0); 5053acd3b1SBarry Smith } 5153acd3b1SBarry Smith 523194b578SJed Brown /* 533194b578SJed Brown Handles setting up the data structure in a call to PetscObjectOptionsBegin() 543194b578SJed Brown */ 554416b707SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionItems *PetscOptionsObject,PetscObject obj) 563194b578SJed Brown { 573194b578SJed Brown PetscErrorCode ierr; 583194b578SJed Brown char title[256]; 593194b578SJed Brown PetscBool flg; 603194b578SJed Brown 613194b578SJed Brown PetscFunctionBegin; 623194b578SJed Brown PetscValidHeader(obj,1); 63e55864a3SBarry Smith PetscOptionsObject->object = obj; 64e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = obj->optionsprinted; 65a297a907SKarl Rupp 663194b578SJed Brown ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr); 673194b578SJed Brown if (flg) { 688caf3d72SBarry Smith ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr); 693194b578SJed Brown } else { 708caf3d72SBarry Smith ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr); 713194b578SJed Brown } 72e55864a3SBarry Smith ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr); 733194b578SJed Brown PetscFunctionReturn(0); 743194b578SJed Brown } 753194b578SJed Brown 7653acd3b1SBarry Smith /* 7753acd3b1SBarry Smith Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd() 7853acd3b1SBarry Smith */ 794416b707SBarry Smith static int PetscOptionItemCreate_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptionItem *amsopt) 8053acd3b1SBarry Smith { 8153acd3b1SBarry Smith int ierr; 824416b707SBarry Smith PetscOptionItem 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 108aee2cecaSBarry Smith /* 1093fc1eb6aSBarry Smith PetscScanString - Gets user input via stdin from process and broadcasts to all processes 1103fc1eb6aSBarry Smith 1113fc1eb6aSBarry Smith Collective on MPI_Comm 1123fc1eb6aSBarry Smith 1133fc1eb6aSBarry Smith Input Parameters: 1143fc1eb6aSBarry Smith + commm - communicator for the broadcast, must be PETSC_COMM_WORLD 1153fc1eb6aSBarry Smith . n - length of the string, must be the same on all processes 1163fc1eb6aSBarry Smith - str - location to store input 117aee2cecaSBarry Smith 118aee2cecaSBarry Smith Bugs: 119aee2cecaSBarry Smith . Assumes process 0 of the given communicator has access to stdin 120aee2cecaSBarry Smith 121aee2cecaSBarry Smith */ 1223fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[]) 123aee2cecaSBarry Smith { 124330cf3c9SBarry Smith size_t i; 125aee2cecaSBarry Smith char c; 1263fc1eb6aSBarry Smith PetscMPIInt rank,nm; 127aee2cecaSBarry Smith PetscErrorCode ierr; 128aee2cecaSBarry Smith 129aee2cecaSBarry Smith PetscFunctionBegin; 130aee2cecaSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 131aee2cecaSBarry Smith if (!rank) { 132aee2cecaSBarry Smith c = (char) getchar(); 133aee2cecaSBarry Smith i = 0; 134aee2cecaSBarry Smith while (c != '\n' && i < n-1) { 135aee2cecaSBarry Smith str[i++] = c; 136aee2cecaSBarry Smith c = (char)getchar(); 137aee2cecaSBarry Smith } 138aee2cecaSBarry Smith str[i] = 0; 139aee2cecaSBarry Smith } 1404dc2109aSBarry Smith ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr); 1413fc1eb6aSBarry Smith ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr); 142aee2cecaSBarry Smith PetscFunctionReturn(0); 143aee2cecaSBarry Smith } 144aee2cecaSBarry Smith 1455b02f95dSBarry Smith /* 1465b02f95dSBarry Smith This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy() 1475b02f95dSBarry Smith */ 1485b02f95dSBarry Smith static PetscErrorCode PetscStrdup(const char s[],char *t[]) 1495b02f95dSBarry Smith { 1505b02f95dSBarry Smith PetscErrorCode ierr; 1515b02f95dSBarry Smith size_t len; 1525b02f95dSBarry Smith char *tmp = 0; 1535b02f95dSBarry Smith 1545b02f95dSBarry Smith PetscFunctionBegin; 1555b02f95dSBarry Smith if (s) { 1565b02f95dSBarry Smith ierr = PetscStrlen(s,&len);CHKERRQ(ierr); 157f416af30SBarry Smith tmp = (char*) malloc((len+1)*sizeof(char)); 1585b02f95dSBarry Smith if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string"); 1595b02f95dSBarry Smith ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr); 1605b02f95dSBarry Smith } 1615b02f95dSBarry Smith *t = tmp; 1625b02f95dSBarry Smith PetscFunctionReturn(0); 1635b02f95dSBarry Smith } 1645b02f95dSBarry Smith 1655b02f95dSBarry Smith 166aee2cecaSBarry Smith /* 1673cc1e11dSBarry Smith PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime 168aee2cecaSBarry Smith 169aee2cecaSBarry Smith Notes: this isn't really practical, it is just to demonstrate the principle 170aee2cecaSBarry Smith 1717781c08eSBarry Smith A carriage return indicates no change from the default; but this like -ksp_monitor <stdout> the default is actually not stdout the default 1727781c08eSBarry Smith is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug? 1737781c08eSBarry Smith 174aee2cecaSBarry Smith Bugs: 1757781c08eSBarry Smith + All processes must traverse through the exact same set of option queries due to the call to PetscScanString() 1763cc1e11dSBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 177aee2cecaSBarry Smith - Only works for PetscInt == int, PetscReal == double etc 178aee2cecaSBarry Smith 1793cc1e11dSBarry Smith Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different 1803cc1e11dSBarry Smith address space and communicating with the PETSc program 1813cc1e11dSBarry Smith 182aee2cecaSBarry Smith */ 1834416b707SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionItems *PetscOptionsObject) 1846356e834SBarry Smith { 1856356e834SBarry Smith PetscErrorCode ierr; 1864416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 1876356e834SBarry Smith char str[512]; 1887781c08eSBarry Smith PetscBool bid; 189a4404d99SBarry Smith PetscReal ir,*valr; 190330cf3c9SBarry Smith PetscInt *vald; 191330cf3c9SBarry Smith size_t i; 1926356e834SBarry Smith 1932a409bb0SBarry Smith PetscFunctionBegin; 194e55864a3SBarry Smith ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr); 1956356e834SBarry Smith while (next) { 1966356e834SBarry Smith switch (next->type) { 1976356e834SBarry Smith case OPTION_HEAD: 1986356e834SBarry Smith break; 199e26ddf31SBarry Smith case OPTION_INT_ARRAY: 200e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 201e26ddf31SBarry Smith vald = (PetscInt*) next->data; 202e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 203e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr); 204e26ddf31SBarry Smith if (i < next->arraylength-1) { 205e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 206e26ddf31SBarry Smith } 207e26ddf31SBarry Smith } 208e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 209e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 210e26ddf31SBarry Smith if (str[0]) { 211e26ddf31SBarry Smith PetscToken token; 212e26ddf31SBarry Smith PetscInt n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end; 213e26ddf31SBarry Smith size_t len; 214e26ddf31SBarry Smith char *value; 215ace3abfcSBarry Smith PetscBool foundrange; 216e26ddf31SBarry Smith 217e26ddf31SBarry Smith next->set = PETSC_TRUE; 218e26ddf31SBarry Smith value = str; 219e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 220e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 221e26ddf31SBarry Smith while (n < nmax) { 222e26ddf31SBarry Smith if (!value) break; 223e26ddf31SBarry Smith 224e26ddf31SBarry Smith /* look for form d-D where d and D are integers */ 225e26ddf31SBarry Smith foundrange = PETSC_FALSE; 226e26ddf31SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 227e26ddf31SBarry Smith if (value[0] == '-') i=2; 228e26ddf31SBarry Smith else i=1; 229330cf3c9SBarry Smith for (;i<len; i++) { 230e26ddf31SBarry Smith if (value[i] == '-') { 231e32f2f54SBarry Smith if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value); 232e26ddf31SBarry Smith value[i] = 0; 233cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr); 234cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr); 235e32f2f54SBarry 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); 236e32f2f54SBarry 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); 237e26ddf31SBarry Smith for (; start<end; start++) { 238e26ddf31SBarry Smith *dvalue = start; dvalue++;n++; 239e26ddf31SBarry Smith } 240e26ddf31SBarry Smith foundrange = PETSC_TRUE; 241e26ddf31SBarry Smith break; 242e26ddf31SBarry Smith } 243e26ddf31SBarry Smith } 244e26ddf31SBarry Smith if (!foundrange) { 245cfbddea1SSatish Balay ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr); 246e26ddf31SBarry Smith dvalue++; 247e26ddf31SBarry Smith n++; 248e26ddf31SBarry Smith } 249e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 250e26ddf31SBarry Smith } 2518c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 252e26ddf31SBarry Smith } 253e26ddf31SBarry Smith break; 254e26ddf31SBarry Smith case OPTION_REAL_ARRAY: 255e55864a3SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr); 256e26ddf31SBarry Smith valr = (PetscReal*) next->data; 257e26ddf31SBarry Smith for (i=0; i<next->arraylength; i++) { 258e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr); 259e26ddf31SBarry Smith if (i < next->arraylength-1) { 260e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr); 261e26ddf31SBarry Smith } 262e26ddf31SBarry Smith } 263e26ddf31SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr); 264e26ddf31SBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 265e26ddf31SBarry Smith if (str[0]) { 266e26ddf31SBarry Smith PetscToken token; 267e26ddf31SBarry Smith PetscInt n = 0,nmax = next->arraylength; 268e26ddf31SBarry Smith PetscReal *dvalue = (PetscReal*)next->data; 269e26ddf31SBarry Smith char *value; 270e26ddf31SBarry Smith 271e26ddf31SBarry Smith next->set = PETSC_TRUE; 272e26ddf31SBarry Smith value = str; 273e26ddf31SBarry Smith ierr = PetscTokenCreate(value,',',&token);CHKERRQ(ierr); 274e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 275e26ddf31SBarry Smith while (n < nmax) { 276e26ddf31SBarry Smith if (!value) break; 277cfbddea1SSatish Balay ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr); 278e26ddf31SBarry Smith dvalue++; 279e26ddf31SBarry Smith n++; 280e26ddf31SBarry Smith ierr = PetscTokenFind(token,&value);CHKERRQ(ierr); 281e26ddf31SBarry Smith } 2828c74ee41SBarry Smith ierr = PetscTokenDestroy(&token);CHKERRQ(ierr); 283e26ddf31SBarry Smith } 284e26ddf31SBarry Smith break; 2856356e834SBarry Smith case OPTION_INT: 286e55864a3SBarry 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); 2873fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 2883fc1eb6aSBarry Smith if (str[0]) { 289d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG) 290d25d7f95SJed Brown long long lid; 291d25d7f95SJed Brown sscanf(str,"%lld",&lid); 2921a1499c8SBarry 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); 293c272547aSJed Brown #else 294d25d7f95SJed Brown long lid; 295d25d7f95SJed Brown sscanf(str,"%ld",&lid); 2961a1499c8SBarry 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); 297c272547aSJed Brown #endif 298a297a907SKarl Rupp 299d25d7f95SJed Brown next->set = PETSC_TRUE; 300d25d7f95SJed Brown *((PetscInt*)next->data) = (PetscInt)lid; 301aee2cecaSBarry Smith } 302aee2cecaSBarry Smith break; 303aee2cecaSBarry Smith case OPTION_REAL: 304e55864a3SBarry 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); 3053fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3063fc1eb6aSBarry Smith if (str[0]) { 307ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 308a4404d99SBarry Smith sscanf(str,"%e",&ir); 309570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 310570b7f6dSBarry Smith float irtemp; 311570b7f6dSBarry Smith sscanf(str,"%e",&irtemp); 312570b7f6dSBarry Smith ir = irtemp; 313ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE) 314aee2cecaSBarry Smith sscanf(str,"%le",&ir); 315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 316d9822059SBarry Smith ir = strtoflt128(str,0); 317d9822059SBarry Smith #else 318513dbe71SLisandro Dalcin SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type"); 319a4404d99SBarry Smith #endif 320aee2cecaSBarry Smith next->set = PETSC_TRUE; 321aee2cecaSBarry Smith *((PetscReal*)next->data) = ir; 322aee2cecaSBarry Smith } 323aee2cecaSBarry Smith break; 3247781c08eSBarry Smith case OPTION_BOOL: 32583355fc5SBarry 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); 3267781c08eSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3277781c08eSBarry Smith if (str[0]) { 3287781c08eSBarry Smith ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr); 3297781c08eSBarry Smith next->set = PETSC_TRUE; 3307781c08eSBarry Smith *((PetscBool*)next->data) = bid; 3317781c08eSBarry Smith } 3327781c08eSBarry Smith break; 333aee2cecaSBarry Smith case OPTION_STRING: 334e55864a3SBarry 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); 3353fc1eb6aSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3363fc1eb6aSBarry Smith if (str[0]) { 337aee2cecaSBarry Smith next->set = PETSC_TRUE; 33864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3395b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3406356e834SBarry Smith } 3416356e834SBarry Smith break; 342a264d7a6SBarry Smith case OPTION_FLIST: 343e55864a3SBarry Smith ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr); 3443cc1e11dSBarry Smith ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr); 3453cc1e11dSBarry Smith if (str[0]) { 346e55864a3SBarry Smith PetscOptionsObject->changedmethod = PETSC_TRUE; 3473cc1e11dSBarry Smith next->set = PETSC_TRUE; 34864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 3495b02f95dSBarry Smith ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr); 3503cc1e11dSBarry Smith } 3513cc1e11dSBarry Smith break; 352b432afdaSMatthew Knepley default: 353b432afdaSMatthew Knepley break; 3546356e834SBarry Smith } 3556356e834SBarry Smith next = next->next; 3566356e834SBarry Smith } 3576356e834SBarry Smith PetscFunctionReturn(0); 3586356e834SBarry Smith } 3596356e834SBarry Smith 360e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 361e04113cfSBarry Smith #include <petscviewersaws.h> 362d5649816SBarry Smith 363d5649816SBarry Smith static int count = 0; 364d5649816SBarry Smith 365e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void) 366d5649816SBarry Smith { 3672657e9d9SBarry Smith PetscFunctionBegin; 368d5649816SBarry Smith PetscFunctionReturn(0); 369d5649816SBarry Smith } 370d5649816SBarry Smith 3719c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n" 37223a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n" 37323a1ff79SBarry Smith "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n" 374d1fc0251SBarry Smith "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n" 37564bbc9efSBarry Smith "<script>\n" 37664bbc9efSBarry Smith "jQuery(document).ready(function() {\n" 37764bbc9efSBarry Smith "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n" 37864bbc9efSBarry Smith "})\n" 37964bbc9efSBarry Smith "</script>\n" 38064bbc9efSBarry Smith "</head>\n"; 3811423471aSBarry Smith 3821423471aSBarry Smith /* Determines the size and style of the scroll region where PETSc options selectable from users are displayed */ 3831423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>"; 38464bbc9efSBarry Smith 385b3506946SBarry Smith /* 3867781c08eSBarry Smith PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs 387b3506946SBarry Smith 388b3506946SBarry Smith Bugs: 389b3506946SBarry Smith + All processes must traverse through the exact same set of option queries do to the call to PetscScanString() 390b3506946SBarry Smith . Internal strings have arbitrary length and string copies are not checked that they fit into string space 391b3506946SBarry Smith - Only works for PetscInt == int, PetscReal == double etc 392b3506946SBarry Smith 393b3506946SBarry Smith 394b3506946SBarry Smith */ 3954416b707SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptionItems *PetscOptionsObject) 396b3506946SBarry Smith { 397b3506946SBarry Smith PetscErrorCode ierr; 3984416b707SBarry Smith PetscOptionItem next = PetscOptionsObject->next; 399d5649816SBarry Smith static int mancount = 0; 400b3506946SBarry Smith char options[16]; 4017aab2a10SBarry Smith PetscBool changedmethod = PETSC_FALSE; 402a23eb982SSurtai Han PetscBool stopasking = PETSC_FALSE; 40388a9752cSBarry Smith char manname[16],textname[16]; 4042657e9d9SBarry Smith char dir[1024]; 405b3506946SBarry Smith 4062a409bb0SBarry Smith PetscFunctionBegin; 407b3506946SBarry Smith /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */ 408b3506946SBarry Smith sprintf(options,"Options_%d",count++); 409a297a907SKarl Rupp 4107325c4a4SBarry Smith PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */ 4111bc75a8dSBarry Smith 412d50831c4SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr); 41383355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING)); 4147781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr); 41583355fc5SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING)); 4162657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN)); 417a23eb982SSurtai Han PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN)); 418b3506946SBarry Smith 419b3506946SBarry Smith while (next) { 420d50831c4SBarry Smith sprintf(manname,"_man_%d",mancount); 4212657e9d9SBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr); 4227781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING)); 423d50831c4SBarry Smith sprintf(textname,"_text_%d",mancount++); 4247781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr); 4257781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING)); 4269f32e415SBarry Smith 427b3506946SBarry Smith switch (next->type) { 428b3506946SBarry Smith case OPTION_HEAD: 429b3506946SBarry Smith break; 430b3506946SBarry Smith case OPTION_INT_ARRAY: 4317781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4322657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT)); 433b3506946SBarry Smith break; 434b3506946SBarry Smith case OPTION_REAL_ARRAY: 4357781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4362657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE)); 437b3506946SBarry Smith break; 438b3506946SBarry Smith case OPTION_INT: 4397781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4402657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT)); 441b3506946SBarry Smith break; 442b3506946SBarry Smith case OPTION_REAL: 4437781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4442657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE)); 445b3506946SBarry Smith break; 4467781c08eSBarry Smith case OPTION_BOOL: 4477781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4482657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN)); 4491ae3d29cSBarry Smith break; 4507781c08eSBarry Smith case OPTION_BOOL_ARRAY: 4517781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4522657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN)); 45371f08665SBarry Smith break; 454b3506946SBarry Smith case OPTION_STRING: 4557781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4567781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 4571ae3d29cSBarry Smith break; 4581ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 4597781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4602657e9d9SBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING)); 461b3506946SBarry Smith break; 462a264d7a6SBarry Smith case OPTION_FLIST: 463a264d7a6SBarry Smith { 464a264d7a6SBarry Smith PetscInt ntext; 4657781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4667781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 467a264d7a6SBarry Smith ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr); 468a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 469a264d7a6SBarry Smith } 470a264d7a6SBarry Smith break; 4711ae3d29cSBarry Smith case OPTION_ELIST: 472a264d7a6SBarry Smith { 473a264d7a6SBarry Smith PetscInt ntext = next->nlist; 4747781c08eSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 4757781c08eSBarry Smith PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING)); 476ead66b60SBarry Smith ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr); 4771ae3d29cSBarry Smith ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr); 478a264d7a6SBarry Smith PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata)); 479a264d7a6SBarry Smith } 480a264d7a6SBarry Smith break; 481b3506946SBarry Smith default: 482b3506946SBarry Smith break; 483b3506946SBarry Smith } 484b3506946SBarry Smith next = next->next; 485b3506946SBarry Smith } 486b3506946SBarry Smith 487b3506946SBarry Smith /* wait until accessor has unlocked the memory */ 48864bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader)); 48964bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom)); 4907aab2a10SBarry Smith ierr = PetscSAWsBlock();CHKERRQ(ierr); 49164bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Header,("index.html")); 49264bbc9efSBarry Smith PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2)); 493b3506946SBarry Smith 49488a9752cSBarry Smith /* determine if any values have been set in GUI */ 49583355fc5SBarry Smith next = PetscOptionsObject->next; 49688a9752cSBarry Smith while (next) { 49788a9752cSBarry Smith ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr); 498f7b25cf6SBarry Smith PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set)); 49988a9752cSBarry Smith next = next->next; 50088a9752cSBarry Smith } 50188a9752cSBarry Smith 502b3506946SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 503f7b25cf6SBarry Smith if (changedmethod) PetscOptionsObject->count = -2; 504b3506946SBarry Smith 505a23eb982SSurtai Han if (stopasking) { 506a23eb982SSurtai Han PetscOptionsPublish = PETSC_FALSE; 50712655325SBarry Smith PetscOptionsObject->count = 0;//do not ask for same thing again 508a23eb982SSurtai Han } 509a23eb982SSurtai Han 5109a492a5cSBarry Smith PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options")); 511b3506946SBarry Smith PetscFunctionReturn(0); 512b3506946SBarry Smith } 513b3506946SBarry Smith #endif 514b3506946SBarry Smith 5154416b707SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionItems *PetscOptionsObject) 51653acd3b1SBarry Smith { 51753acd3b1SBarry Smith PetscErrorCode ierr; 5184416b707SBarry Smith PetscOptionItem last; 5196356e834SBarry Smith char option[256],value[1024],tmp[32]; 520330cf3c9SBarry Smith size_t j; 52153acd3b1SBarry Smith 52253acd3b1SBarry Smith PetscFunctionBegin; 52383355fc5SBarry Smith if (PetscOptionsObject->next) { 52483355fc5SBarry Smith if (!PetscOptionsObject->count) { 525a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS) 526f7b25cf6SBarry Smith ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr); 527b3506946SBarry Smith #else 528e55864a3SBarry Smith ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr); 529b3506946SBarry Smith #endif 530aee2cecaSBarry Smith } 531aee2cecaSBarry Smith } 5326356e834SBarry Smith 533e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr); 5346356e834SBarry Smith 535e26ddf31SBarry Smith /* reset counter to -2; this updates the screen with the new options for the selected method */ 536e55864a3SBarry Smith if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2; 5377a72a596SBarry Smith /* reset alreadyprinted flag */ 538e55864a3SBarry Smith PetscOptionsObject->alreadyprinted = PETSC_FALSE; 539e55864a3SBarry Smith if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE; 540e55864a3SBarry Smith PetscOptionsObject->object = NULL; 54153acd3b1SBarry Smith 542e55864a3SBarry Smith while (PetscOptionsObject->next) { 543e55864a3SBarry Smith if (PetscOptionsObject->next->set) { 544e55864a3SBarry Smith if (PetscOptionsObject->prefix) { 54553acd3b1SBarry Smith ierr = PetscStrcpy(option,"-");CHKERRQ(ierr); 546e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr); 547e55864a3SBarry Smith ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr); 5486356e834SBarry Smith } else { 549e55864a3SBarry Smith ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr); 5506356e834SBarry Smith } 5516356e834SBarry Smith 552e55864a3SBarry Smith switch (PetscOptionsObject->next->type) { 5536356e834SBarry Smith case OPTION_HEAD: 5546356e834SBarry Smith break; 5556356e834SBarry Smith case OPTION_INT_ARRAY: 556e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]); 557e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 558e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]); 5596356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5606356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5616356e834SBarry Smith } 5626356e834SBarry Smith break; 5636356e834SBarry Smith case OPTION_INT: 564e55864a3SBarry Smith sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data); 5656356e834SBarry Smith break; 5666356e834SBarry Smith case OPTION_REAL: 567e55864a3SBarry Smith sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data); 5686356e834SBarry Smith break; 5696356e834SBarry Smith case OPTION_REAL_ARRAY: 570e55864a3SBarry Smith sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]); 571e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 572e55864a3SBarry Smith sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]); 5736356e834SBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5746356e834SBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5756356e834SBarry Smith } 5766356e834SBarry Smith break; 577050cccc3SHong Zhang case OPTION_SCALAR_ARRAY: 57895f3a755SJose E. Roman sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[0])); 579050cccc3SHong Zhang for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 58095f3a755SJose E. Roman sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscScalar*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscScalar*)PetscOptionsObject->next->data)[j])); 581050cccc3SHong Zhang ierr = PetscStrcat(value,",");CHKERRQ(ierr); 582050cccc3SHong Zhang ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 583050cccc3SHong Zhang } 584050cccc3SHong Zhang break; 5857781c08eSBarry Smith case OPTION_BOOL: 586e55864a3SBarry Smith sprintf(value,"%d",*(int*)PetscOptionsObject->next->data); 5876356e834SBarry Smith break; 5887781c08eSBarry Smith case OPTION_BOOL_ARRAY: 589e55864a3SBarry Smith sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]); 590e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 591e55864a3SBarry Smith sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]); 5921ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 5931ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 5941ae3d29cSBarry Smith } 5951ae3d29cSBarry Smith break; 596a264d7a6SBarry Smith case OPTION_FLIST: 5976991f827SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 5986991f827SBarry Smith break; 5991ae3d29cSBarry Smith case OPTION_ELIST: 600e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 6016356e834SBarry Smith break; 6021ae3d29cSBarry Smith case OPTION_STRING: 603e55864a3SBarry Smith ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr); 604d51da6bfSBarry Smith break; 6051ae3d29cSBarry Smith case OPTION_STRING_ARRAY: 606e55864a3SBarry Smith sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]); 607e55864a3SBarry Smith for (j=1; j<PetscOptionsObject->next->arraylength; j++) { 608e55864a3SBarry Smith sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]); 6091ae3d29cSBarry Smith ierr = PetscStrcat(value,",");CHKERRQ(ierr); 6101ae3d29cSBarry Smith ierr = PetscStrcat(value,tmp);CHKERRQ(ierr); 6111ae3d29cSBarry Smith } 6126356e834SBarry Smith break; 6136356e834SBarry Smith } 614c5929fdfSBarry Smith ierr = PetscOptionsSetValue(PetscOptionsObject->options,option,value);CHKERRQ(ierr); 6156356e834SBarry Smith } 6166991f827SBarry Smith if (PetscOptionsObject->next->type == OPTION_ELIST) { 6176991f827SBarry Smith ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr); 6186991f827SBarry Smith } 619e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr); 620e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr); 621e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr); 622e55864a3SBarry Smith ierr = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr); 623c979a496SBarry Smith 62483355fc5SBarry Smith if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){ 62583355fc5SBarry Smith free(PetscOptionsObject->next->data); 626c979a496SBarry Smith } else { 62783355fc5SBarry Smith ierr = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr); 628c979a496SBarry Smith } 6297781c08eSBarry Smith 63083355fc5SBarry Smith last = PetscOptionsObject->next; 63183355fc5SBarry Smith PetscOptionsObject->next = PetscOptionsObject->next->next; 6326356e834SBarry Smith ierr = PetscFree(last);CHKERRQ(ierr); 6336356e834SBarry Smith } 634f59f755dSBarry Smith ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr); 635e55864a3SBarry Smith PetscOptionsObject->next = 0; 63653acd3b1SBarry Smith PetscFunctionReturn(0); 63753acd3b1SBarry Smith } 63853acd3b1SBarry Smith 63953acd3b1SBarry Smith /*@C 64053acd3b1SBarry Smith PetscOptionsEnum - Gets the enum value for a particular option in the database. 64153acd3b1SBarry Smith 6423f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 64353acd3b1SBarry Smith 64453acd3b1SBarry Smith Input Parameters: 64553acd3b1SBarry Smith + opt - option name 64653acd3b1SBarry Smith . text - short string that describes the option 64753acd3b1SBarry Smith . man - manual page with additional information on option 64853acd3b1SBarry Smith . list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null 6490fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 6500fdccdaeSBarry Smith $ PetscOptionsEnum(..., obj->value,&object->value,...) or 6510fdccdaeSBarry Smith $ value = defaultvalue 6520fdccdaeSBarry Smith $ PetscOptionsEnum(..., value,&value,&flg); 6530fdccdaeSBarry Smith $ if (flg) { 65453acd3b1SBarry Smith 65553acd3b1SBarry Smith Output Parameter: 65653acd3b1SBarry Smith + value - the value to return 657b32e0204SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 65853acd3b1SBarry Smith 65953acd3b1SBarry Smith Level: beginner 66053acd3b1SBarry Smith 66153acd3b1SBarry Smith Concepts: options database 66253acd3b1SBarry Smith 66353acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 66453acd3b1SBarry Smith 66553acd3b1SBarry Smith list is usually something like PCASMTypes or some other predefined list of enum names 66653acd3b1SBarry Smith 667*2efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 668*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 669*2efd9cb1SBarry Smith 67053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 671acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 672acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 67353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 67453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 675acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 676a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 67753acd3b1SBarry Smith @*/ 6784416b707SBarry Smith PetscErrorCode PetscOptionsEnum_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool *set) 67953acd3b1SBarry Smith { 68053acd3b1SBarry Smith PetscErrorCode ierr; 68153acd3b1SBarry Smith PetscInt ntext = 0; 682aa5bb8c0SSatish Balay PetscInt tval; 683ace3abfcSBarry Smith PetscBool tflg; 68453acd3b1SBarry Smith 68553acd3b1SBarry Smith PetscFunctionBegin; 68653acd3b1SBarry Smith while (list[ntext++]) { 687e32f2f54SBarry Smith if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries"); 68853acd3b1SBarry Smith } 689e32f2f54SBarry Smith if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 69053acd3b1SBarry Smith ntext -= 3; 691e55864a3SBarry Smith ierr = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr); 692aa5bb8c0SSatish Balay /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */ 693aa5bb8c0SSatish Balay if (tflg) *value = (PetscEnum)tval; 694aa5bb8c0SSatish Balay if (set) *set = tflg; 69553acd3b1SBarry Smith PetscFunctionReturn(0); 69653acd3b1SBarry Smith } 69753acd3b1SBarry Smith 698d3e47460SLisandro Dalcin /*@C 699d3e47460SLisandro Dalcin PetscOptionsEnumArray - Gets an array of enum values for a particular 700d3e47460SLisandro Dalcin option in the database. 701d3e47460SLisandro Dalcin 702d3e47460SLisandro Dalcin Logically Collective on the communicator passed in PetscOptionsBegin() 703d3e47460SLisandro Dalcin 704d3e47460SLisandro Dalcin Input Parameters: 705d3e47460SLisandro Dalcin + opt - the option one is seeking 706d3e47460SLisandro Dalcin . text - short string describing option 707d3e47460SLisandro Dalcin . man - manual page for option 708d3e47460SLisandro Dalcin - n - maximum number of values 709d3e47460SLisandro Dalcin 710d3e47460SLisandro Dalcin Output Parameter: 711d3e47460SLisandro Dalcin + value - location to copy values 712d3e47460SLisandro Dalcin . n - actual number of values found 713d3e47460SLisandro Dalcin - set - PETSC_TRUE if found, else PETSC_FALSE 714d3e47460SLisandro Dalcin 715d3e47460SLisandro Dalcin Level: beginner 716d3e47460SLisandro Dalcin 717d3e47460SLisandro Dalcin Notes: 718d3e47460SLisandro Dalcin The array must be passed as a comma separated list. 719d3e47460SLisandro Dalcin 720d3e47460SLisandro Dalcin There must be no intervening spaces between the values. 721d3e47460SLisandro Dalcin 722d3e47460SLisandro Dalcin Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 723d3e47460SLisandro Dalcin 724d3e47460SLisandro Dalcin Concepts: options database^array of enums 725d3e47460SLisandro Dalcin 726d3e47460SLisandro Dalcin .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(), 727d3e47460SLisandro Dalcin PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 728d3e47460SLisandro Dalcin PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 729d3e47460SLisandro Dalcin PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 730d3e47460SLisandro Dalcin PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 731d3e47460SLisandro Dalcin PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray() 732d3e47460SLisandro Dalcin @*/ 7334416b707SBarry Smith PetscErrorCode PetscOptionsEnumArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum value[],PetscInt *n,PetscBool *set) 734d3e47460SLisandro Dalcin { 735d3e47460SLisandro Dalcin PetscInt i,nlist = 0; 7364416b707SBarry Smith PetscOptionItem amsopt; 737d3e47460SLisandro Dalcin PetscErrorCode ierr; 738d3e47460SLisandro Dalcin 739d3e47460SLisandro Dalcin PetscFunctionBegin; 740d3e47460SLisandro 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"); 741d3e47460SLisandro Dalcin if (nlist < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix"); 742d3e47460SLisandro Dalcin nlist -= 3; /* drop enum name, prefix, and null termination */ 743d3e47460SLisandro Dalcin if (0 && !PetscOptionsObject->count) { /* XXX Requires additional support */ 744d3e47460SLisandro Dalcin PetscEnum *vals; 7454416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY/*XXX OPTION_ENUM_ARRAY*/,&amsopt);CHKERRQ(ierr); 746d3e47460SLisandro Dalcin ierr = PetscStrNArrayallocpy(nlist,list,(char***)&amsopt->list);CHKERRQ(ierr); 747d3e47460SLisandro Dalcin amsopt->nlist = nlist; 748d3e47460SLisandro Dalcin ierr = PetscMalloc1(*n,(PetscEnum**)&amsopt->data);CHKERRQ(ierr); 749d3e47460SLisandro Dalcin amsopt->arraylength = *n; 750d3e47460SLisandro Dalcin vals = (PetscEnum*)amsopt->data; 751d3e47460SLisandro Dalcin for (i=0; i<*n; i++) vals[i] = value[i]; 752d3e47460SLisandro Dalcin } 753c5929fdfSBarry Smith ierr = PetscOptionsGetEnumArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,value,n,set);CHKERRQ(ierr); 754d3e47460SLisandro Dalcin if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 755d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,list[value[0]]);CHKERRQ(ierr); 756d3e47460SLisandro Dalcin for (i=1; i<*n; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%s",list[value[i]]);CHKERRQ(ierr);} 757d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (choose from)",text);CHKERRQ(ierr); 758d3e47460SLisandro Dalcin for (i=0; i<nlist; i++) {ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);} 759d3e47460SLisandro Dalcin ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr); 760d3e47460SLisandro Dalcin } 761d3e47460SLisandro Dalcin PetscFunctionReturn(0); 762d3e47460SLisandro Dalcin } 763d3e47460SLisandro Dalcin 76453acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/ 76553acd3b1SBarry Smith /*@C 76653acd3b1SBarry Smith PetscOptionsInt - Gets the integer value for a particular option in the database. 76753acd3b1SBarry Smith 7683f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 76953acd3b1SBarry Smith 77053acd3b1SBarry Smith Input Parameters: 77153acd3b1SBarry Smith + opt - option name 77253acd3b1SBarry Smith . text - short string that describes the option 77353acd3b1SBarry Smith . man - manual page with additional information on option 7740fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 7750fdccdaeSBarry Smith $ PetscOptionsInt(..., obj->value,&object->value,...) or 7760fdccdaeSBarry Smith $ value = defaultvalue 7770fdccdaeSBarry Smith $ PetscOptionsInt(..., value,&value,&flg); 7780fdccdaeSBarry Smith $ if (flg) { 77953acd3b1SBarry Smith 78053acd3b1SBarry Smith Output Parameter: 78153acd3b1SBarry Smith + value - the integer value to return 78253acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 78353acd3b1SBarry Smith 784*2efd9cb1SBarry Smith Notes: If the user does not supply the option at all value is NOT changed. Thus 785*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 786*2efd9cb1SBarry Smith 78753acd3b1SBarry Smith Level: beginner 78853acd3b1SBarry Smith 78953acd3b1SBarry Smith Concepts: options database^has int 79053acd3b1SBarry Smith 79153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 79253acd3b1SBarry Smith 79353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(), 794acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 795acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 79653acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 79753acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 798acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 799a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 80053acd3b1SBarry Smith @*/ 8014416b707SBarry Smith PetscErrorCode PetscOptionsInt_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool *set) 80253acd3b1SBarry Smith { 80353acd3b1SBarry Smith PetscErrorCode ierr; 8044416b707SBarry Smith PetscOptionItem amsopt; 80512655325SBarry Smith PetscBool wasset; 80653acd3b1SBarry Smith 80753acd3b1SBarry Smith PetscFunctionBegin; 808e55864a3SBarry Smith if (!PetscOptionsObject->count) { 8094416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr); 8106356e834SBarry Smith ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr); 81112655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 8123e211508SBarry Smith 813c5929fdfSBarry Smith ierr = PetscOptionsGetInt(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,¤tvalue,&wasset);CHKERRQ(ierr); 8143e211508SBarry Smith if (wasset) { 81512655325SBarry Smith *(PetscInt*)amsopt->data = currentvalue; 8163e211508SBarry Smith } 817af6d86caSBarry Smith } 818c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 819e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 8201a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 82153acd3b1SBarry Smith } 82253acd3b1SBarry Smith PetscFunctionReturn(0); 82353acd3b1SBarry Smith } 82453acd3b1SBarry Smith 82553acd3b1SBarry Smith /*@C 82653acd3b1SBarry Smith PetscOptionsString - Gets the string value for a particular option in the database. 82753acd3b1SBarry Smith 8283f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 82953acd3b1SBarry Smith 83053acd3b1SBarry Smith Input Parameters: 83153acd3b1SBarry Smith + opt - option name 83253acd3b1SBarry Smith . text - short string that describes the option 83353acd3b1SBarry Smith . man - manual page with additional information on option 8340fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value 835bcbf2dc5SJed Brown - len - length of the result string including null terminator 83653acd3b1SBarry Smith 83753acd3b1SBarry Smith Output Parameter: 83853acd3b1SBarry Smith + value - the value to return 83953acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 84053acd3b1SBarry Smith 84153acd3b1SBarry Smith Level: beginner 84253acd3b1SBarry Smith 84353acd3b1SBarry Smith Concepts: options database^has int 84453acd3b1SBarry Smith 84553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 84653acd3b1SBarry Smith 8477fccdfe4SBarry 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). 8487fccdfe4SBarry Smith 849*2efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 850*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 851*2efd9cb1SBarry Smith 852*2efd9cb1SBarry Smith 853c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 854acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 855acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(), 85653acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 85753acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 858acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 859a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 86053acd3b1SBarry Smith @*/ 8614416b707SBarry Smith PetscErrorCode PetscOptionsString_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool *set) 86253acd3b1SBarry Smith { 86353acd3b1SBarry Smith PetscErrorCode ierr; 8644416b707SBarry Smith PetscOptionItem amsopt; 86553acd3b1SBarry Smith 86653acd3b1SBarry Smith PetscFunctionBegin; 8671a1499c8SBarry Smith if (!PetscOptionsObject->count) { 8684416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 86964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 8700fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 871af6d86caSBarry Smith } 872c5929fdfSBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 873e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 8741a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr); 87553acd3b1SBarry Smith } 87653acd3b1SBarry Smith PetscFunctionReturn(0); 87753acd3b1SBarry Smith } 87853acd3b1SBarry Smith 87953acd3b1SBarry Smith /*@C 88053acd3b1SBarry Smith PetscOptionsReal - Gets the PetscReal value for a particular option in the database. 88153acd3b1SBarry Smith 8823f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 88353acd3b1SBarry Smith 88453acd3b1SBarry Smith Input Parameters: 88553acd3b1SBarry Smith + opt - option name 88653acd3b1SBarry Smith . text - short string that describes the option 88753acd3b1SBarry Smith . man - manual page with additional information on option 8880fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 8890fdccdaeSBarry Smith $ PetscOptionsReal(..., obj->value,&object->value,...) or 8900fdccdaeSBarry Smith $ value = defaultvalue 8910fdccdaeSBarry Smith $ PetscOptionsReal(..., value,&value,&flg); 8920fdccdaeSBarry Smith $ if (flg) { 89353acd3b1SBarry Smith 89453acd3b1SBarry Smith Output Parameter: 89553acd3b1SBarry Smith + value - the value to return 89653acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 89753acd3b1SBarry Smith 898*2efd9cb1SBarry Smith Notes: If the user does not supply the option at all value is NOT changed. Thus 899*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 900*2efd9cb1SBarry Smith 90153acd3b1SBarry Smith Level: beginner 90253acd3b1SBarry Smith 90353acd3b1SBarry Smith Concepts: options database^has int 90453acd3b1SBarry Smith 90553acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 90653acd3b1SBarry Smith 907c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 908acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 909acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 91053acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 91153acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 912acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 913a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 91453acd3b1SBarry Smith @*/ 9154416b707SBarry Smith PetscErrorCode PetscOptionsReal_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool *set) 91653acd3b1SBarry Smith { 91753acd3b1SBarry Smith PetscErrorCode ierr; 9184416b707SBarry Smith PetscOptionItem amsopt; 91953acd3b1SBarry Smith 92053acd3b1SBarry Smith PetscFunctionBegin; 921e55864a3SBarry Smith if (!PetscOptionsObject->count) { 9224416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr); 923538aa990SBarry Smith ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 924a297a907SKarl Rupp 9250fdccdaeSBarry Smith *(PetscReal*)amsopt->data = currentvalue; 926538aa990SBarry Smith } 927c5929fdfSBarry Smith ierr = PetscOptionsGetReal(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 9281a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 9291a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr); 93053acd3b1SBarry Smith } 93153acd3b1SBarry Smith PetscFunctionReturn(0); 93253acd3b1SBarry Smith } 93353acd3b1SBarry Smith 93453acd3b1SBarry Smith /*@C 93553acd3b1SBarry Smith PetscOptionsScalar - Gets the scalar value for a particular option in the database. 93653acd3b1SBarry Smith 9373f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 93853acd3b1SBarry Smith 93953acd3b1SBarry Smith Input Parameters: 94053acd3b1SBarry Smith + opt - option name 94153acd3b1SBarry Smith . text - short string that describes the option 94253acd3b1SBarry Smith . man - manual page with additional information on option 9430fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either 9440fdccdaeSBarry Smith $ PetscOptionsScalar(..., obj->value,&object->value,...) or 9450fdccdaeSBarry Smith $ value = defaultvalue 9460fdccdaeSBarry Smith $ PetscOptionsScalar(..., value,&value,&flg); 9470fdccdaeSBarry Smith $ if (flg) { 9480fdccdaeSBarry Smith 94953acd3b1SBarry Smith 95053acd3b1SBarry Smith Output Parameter: 95153acd3b1SBarry Smith + value - the value to return 95253acd3b1SBarry Smith - flg - PETSC_TRUE if found, else PETSC_FALSE 95353acd3b1SBarry Smith 954*2efd9cb1SBarry Smith Notes: If the user does not supply the option at all value is NOT changed. Thus 955*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 956*2efd9cb1SBarry Smith 95753acd3b1SBarry Smith Level: beginner 95853acd3b1SBarry Smith 95953acd3b1SBarry Smith Concepts: options database^has int 96053acd3b1SBarry Smith 96153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 96253acd3b1SBarry Smith 963c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 964acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 965acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 96653acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 96753acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 968acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 969a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 97053acd3b1SBarry Smith @*/ 9714416b707SBarry Smith PetscErrorCode PetscOptionsScalar_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool *set) 97253acd3b1SBarry Smith { 97353acd3b1SBarry Smith PetscErrorCode ierr; 97453acd3b1SBarry Smith 97553acd3b1SBarry Smith PetscFunctionBegin; 97653acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 9770fdccdaeSBarry Smith ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr); 97853acd3b1SBarry Smith #else 979c5929fdfSBarry Smith ierr = PetscOptionsGetScalar(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr); 98053acd3b1SBarry Smith #endif 98153acd3b1SBarry Smith PetscFunctionReturn(0); 98253acd3b1SBarry Smith } 98353acd3b1SBarry Smith 98453acd3b1SBarry Smith /*@C 98590d69ab7SBarry 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 98690d69ab7SBarry Smith its value is set to false. 98753acd3b1SBarry Smith 9883f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 98953acd3b1SBarry Smith 99053acd3b1SBarry Smith Input Parameters: 99153acd3b1SBarry Smith + opt - option name 99253acd3b1SBarry Smith . text - short string that describes the option 99353acd3b1SBarry Smith - man - manual page with additional information on option 99453acd3b1SBarry Smith 99553acd3b1SBarry Smith Output Parameter: 99653acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 99753acd3b1SBarry Smith 99853acd3b1SBarry Smith Level: beginner 99953acd3b1SBarry Smith 100053acd3b1SBarry Smith Concepts: options database^has int 100153acd3b1SBarry Smith 100253acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 100353acd3b1SBarry Smith 1004c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 1005acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1006acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 100753acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 100853acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1009acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1010a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 101153acd3b1SBarry Smith @*/ 10124416b707SBarry Smith PetscErrorCode PetscOptionsName_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 101353acd3b1SBarry Smith { 101453acd3b1SBarry Smith PetscErrorCode ierr; 10154416b707SBarry Smith PetscOptionItem amsopt; 101653acd3b1SBarry Smith 101753acd3b1SBarry Smith PetscFunctionBegin; 1018e55864a3SBarry Smith if (!PetscOptionsObject->count) { 10194416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1020ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1021a297a907SKarl Rupp 1022ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 10231ae3d29cSBarry Smith } 1024c5929fdfSBarry Smith ierr = PetscOptionsHasName(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr); 1025e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1026e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 102753acd3b1SBarry Smith } 102853acd3b1SBarry Smith PetscFunctionReturn(0); 102953acd3b1SBarry Smith } 103053acd3b1SBarry Smith 103153acd3b1SBarry Smith /*@C 1032a264d7a6SBarry Smith PetscOptionsFList - Puts a list of option values that a single one may be selected from 103353acd3b1SBarry Smith 10343f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 103553acd3b1SBarry Smith 103653acd3b1SBarry Smith Input Parameters: 103753acd3b1SBarry Smith + opt - option name 103853acd3b1SBarry Smith . text - short string that describes the option 103953acd3b1SBarry Smith . man - manual page with additional information on option 104053acd3b1SBarry Smith . list - the possible choices 10410fdccdaeSBarry Smith . currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 10420fdccdaeSBarry Smith $ PetscOptionsFlist(..., obj->value,value,len,&flg); 10430fdccdaeSBarry Smith $ if (flg) { 10443cc1e11dSBarry Smith - len - the length of the character array value 104553acd3b1SBarry Smith 104653acd3b1SBarry Smith Output Parameter: 104753acd3b1SBarry Smith + value - the value to return 104853acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 104953acd3b1SBarry Smith 105053acd3b1SBarry Smith Level: intermediate 105153acd3b1SBarry Smith 105253acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 105353acd3b1SBarry Smith 1054*2efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 1055*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 1056*2efd9cb1SBarry Smith 105753acd3b1SBarry Smith See PetscOptionsEList() for when the choices are given in a string array 105853acd3b1SBarry Smith 105953acd3b1SBarry Smith To get a listing of all currently specified options, 106088c29154SBarry Smith see PetscOptionsView() or PetscOptionsGetAll() 106153acd3b1SBarry Smith 1062eabe10d7SBarry Smith Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list 1063eabe10d7SBarry Smith 106453acd3b1SBarry Smith Concepts: options database^list 106553acd3b1SBarry Smith 1066c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1067acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 106853acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 106953acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1070acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1071a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum() 107253acd3b1SBarry Smith @*/ 10734416b707SBarry Smith PetscErrorCode PetscOptionsFList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool *set) 107453acd3b1SBarry Smith { 107553acd3b1SBarry Smith PetscErrorCode ierr; 10764416b707SBarry Smith PetscOptionItem amsopt; 107753acd3b1SBarry Smith 107853acd3b1SBarry Smith PetscFunctionBegin; 10791a1499c8SBarry Smith if (!PetscOptionsObject->count) { 10804416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr); 108164facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 10820fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 10833cc1e11dSBarry Smith amsopt->flist = list; 10843cc1e11dSBarry Smith } 1085c5929fdfSBarry Smith ierr = PetscOptionsGetString(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr); 10861a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 10871a1499c8SBarry Smith ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr); 108853acd3b1SBarry Smith } 108953acd3b1SBarry Smith PetscFunctionReturn(0); 109053acd3b1SBarry Smith } 109153acd3b1SBarry Smith 109253acd3b1SBarry Smith /*@C 109353acd3b1SBarry Smith PetscOptionsEList - Puts a list of option values that a single one may be selected from 109453acd3b1SBarry Smith 10953f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 109653acd3b1SBarry Smith 109753acd3b1SBarry Smith Input Parameters: 109853acd3b1SBarry Smith + opt - option name 109953acd3b1SBarry Smith . ltext - short string that describes the option 110053acd3b1SBarry Smith . man - manual page with additional information on option 1101a264d7a6SBarry Smith . list - the possible choices (one of these must be selected, anything else is invalid) 110253acd3b1SBarry Smith . ntext - number of choices 11030fdccdaeSBarry Smith - currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with 11040fdccdaeSBarry Smith $ PetscOptionsElist(..., obj->value,&value,&flg); 11050fdccdaeSBarry Smith $ if (flg) { 11060fdccdaeSBarry Smith 110753acd3b1SBarry Smith 110853acd3b1SBarry Smith Output Parameter: 110953acd3b1SBarry Smith + value - the index of the value to return 111053acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 111153acd3b1SBarry Smith 111253acd3b1SBarry Smith Level: intermediate 111353acd3b1SBarry Smith 111453acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 111553acd3b1SBarry Smith 1116*2efd9cb1SBarry Smith If the user does not supply the option at all value is NOT changed. Thus 1117*2efd9cb1SBarry Smith you should ALWAYS initialize value if you access it without first checking if the set flag is true. 1118*2efd9cb1SBarry Smith 1119a264d7a6SBarry Smith See PetscOptionsFList() for when the choices are given in a PetscFunctionList() 112053acd3b1SBarry Smith 112153acd3b1SBarry Smith Concepts: options database^list 112253acd3b1SBarry Smith 1123c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1124acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 112553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 112653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1127acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1128a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEnum() 112953acd3b1SBarry Smith @*/ 11304416b707SBarry Smith PetscErrorCode PetscOptionsEList_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool *set) 113153acd3b1SBarry Smith { 113253acd3b1SBarry Smith PetscErrorCode ierr; 113353acd3b1SBarry Smith PetscInt i; 11344416b707SBarry Smith PetscOptionItem amsopt; 113553acd3b1SBarry Smith 113653acd3b1SBarry Smith PetscFunctionBegin; 11371a1499c8SBarry Smith if (!PetscOptionsObject->count) { 11384416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr); 113964facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 11400fdccdaeSBarry Smith ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr); 11416991f827SBarry Smith ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr); 11421ae3d29cSBarry Smith amsopt->nlist = ntext; 11431ae3d29cSBarry Smith } 1144c5929fdfSBarry Smith ierr = PetscOptionsGetEList(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr); 11451a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1146e3f729a5SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s> %s (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue,ltext);CHKERRQ(ierr); 114753acd3b1SBarry Smith for (i=0; i<ntext; i++) { 1148e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr); 114953acd3b1SBarry Smith } 1150e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr); 115153acd3b1SBarry Smith } 115253acd3b1SBarry Smith PetscFunctionReturn(0); 115353acd3b1SBarry Smith } 115453acd3b1SBarry Smith 115553acd3b1SBarry Smith /*@C 1156acfcf0e5SJed Brown PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for 1157d5649816SBarry Smith which at most a single value can be true. 115853acd3b1SBarry Smith 11593f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 116053acd3b1SBarry Smith 116153acd3b1SBarry Smith Input Parameters: 116253acd3b1SBarry Smith + opt - option name 116353acd3b1SBarry Smith . text - short string that describes the option 116453acd3b1SBarry Smith - man - manual page with additional information on option 116553acd3b1SBarry Smith 116653acd3b1SBarry Smith Output Parameter: 116753acd3b1SBarry Smith . flg - whether that option was set or not 116853acd3b1SBarry Smith 116953acd3b1SBarry Smith Level: intermediate 117053acd3b1SBarry Smith 117153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 117253acd3b1SBarry Smith 1173acfcf0e5SJed Brown Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd() 117453acd3b1SBarry Smith 117553acd3b1SBarry Smith Concepts: options database^logical group 117653acd3b1SBarry Smith 1177c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1178acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 117953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 118053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1181acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1182a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 118353acd3b1SBarry Smith @*/ 11844416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroupBegin_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 118553acd3b1SBarry Smith { 118653acd3b1SBarry Smith PetscErrorCode ierr; 11874416b707SBarry Smith PetscOptionItem amsopt; 118853acd3b1SBarry Smith 118953acd3b1SBarry Smith PetscFunctionBegin; 1190e55864a3SBarry Smith if (!PetscOptionsObject->count) { 11914416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1192ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1193a297a907SKarl Rupp 1194ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 11951ae3d29cSBarry Smith } 119668b16fdaSBarry Smith *flg = PETSC_FALSE; 1197c5929fdfSBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1198e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1199e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," Pick at most one of -------------\n");CHKERRQ(ierr); 1200e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 120153acd3b1SBarry Smith } 120253acd3b1SBarry Smith PetscFunctionReturn(0); 120353acd3b1SBarry Smith } 120453acd3b1SBarry Smith 120553acd3b1SBarry Smith /*@C 1206acfcf0e5SJed Brown PetscOptionsBoolGroup - One in a series of logical queries on the options database for 1207d5649816SBarry Smith which at most a single value can be true. 120853acd3b1SBarry Smith 12093f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 121053acd3b1SBarry Smith 121153acd3b1SBarry Smith Input Parameters: 121253acd3b1SBarry Smith + opt - option name 121353acd3b1SBarry Smith . text - short string that describes the option 121453acd3b1SBarry Smith - man - manual page with additional information on option 121553acd3b1SBarry Smith 121653acd3b1SBarry Smith Output Parameter: 121753acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 121853acd3b1SBarry Smith 121953acd3b1SBarry Smith Level: intermediate 122053acd3b1SBarry Smith 122153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 122253acd3b1SBarry Smith 1223acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd() 122453acd3b1SBarry Smith 122553acd3b1SBarry Smith Concepts: options database^logical group 122653acd3b1SBarry Smith 1227c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1228acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 122953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 123053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1231acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1232a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 123353acd3b1SBarry Smith @*/ 12344416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroup_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 123553acd3b1SBarry Smith { 123653acd3b1SBarry Smith PetscErrorCode ierr; 12374416b707SBarry Smith PetscOptionItem amsopt; 123853acd3b1SBarry Smith 123953acd3b1SBarry Smith PetscFunctionBegin; 1240e55864a3SBarry Smith if (!PetscOptionsObject->count) { 12414416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1242ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1243a297a907SKarl Rupp 1244ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 12451ae3d29cSBarry Smith } 124617326d04SJed Brown *flg = PETSC_FALSE; 1247c5929fdfSBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1248e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1249e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 125053acd3b1SBarry Smith } 125153acd3b1SBarry Smith PetscFunctionReturn(0); 125253acd3b1SBarry Smith } 125353acd3b1SBarry Smith 125453acd3b1SBarry Smith /*@C 1255acfcf0e5SJed Brown PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for 1256d5649816SBarry Smith which at most a single value can be true. 125753acd3b1SBarry Smith 12583f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 125953acd3b1SBarry Smith 126053acd3b1SBarry Smith Input Parameters: 126153acd3b1SBarry Smith + opt - option name 126253acd3b1SBarry Smith . text - short string that describes the option 126353acd3b1SBarry Smith - man - manual page with additional information on option 126453acd3b1SBarry Smith 126553acd3b1SBarry Smith Output Parameter: 126653acd3b1SBarry Smith . flg - PETSC_TRUE if found, else PETSC_FALSE 126753acd3b1SBarry Smith 126853acd3b1SBarry Smith Level: intermediate 126953acd3b1SBarry Smith 127053acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 127153acd3b1SBarry Smith 1272acfcf0e5SJed Brown Must follow a PetscOptionsBoolGroupBegin() 127353acd3b1SBarry Smith 127453acd3b1SBarry Smith Concepts: options database^logical group 127553acd3b1SBarry Smith 1276c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1277acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 127853acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 127953acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1280acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1281a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 128253acd3b1SBarry Smith @*/ 12834416b707SBarry Smith PetscErrorCode PetscOptionsBoolGroupEnd_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool *flg) 128453acd3b1SBarry Smith { 128553acd3b1SBarry Smith PetscErrorCode ierr; 12864416b707SBarry Smith PetscOptionItem amsopt; 128753acd3b1SBarry Smith 128853acd3b1SBarry Smith PetscFunctionBegin; 1289e55864a3SBarry Smith if (!PetscOptionsObject->count) { 12904416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1291ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1292a297a907SKarl Rupp 1293ace3abfcSBarry Smith *(PetscBool*)amsopt->data = PETSC_FALSE; 12941ae3d29cSBarry Smith } 129517326d04SJed Brown *flg = PETSC_FALSE; 1296c5929fdfSBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr); 1297e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1298e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 129953acd3b1SBarry Smith } 130053acd3b1SBarry Smith PetscFunctionReturn(0); 130153acd3b1SBarry Smith } 130253acd3b1SBarry Smith 130353acd3b1SBarry Smith /*@C 1304acfcf0e5SJed Brown PetscOptionsBool - Determines if a particular option is in the database with a true or false 130553acd3b1SBarry Smith 13063f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 130753acd3b1SBarry Smith 130853acd3b1SBarry Smith Input Parameters: 130953acd3b1SBarry Smith + opt - option name 131053acd3b1SBarry Smith . text - short string that describes the option 1311868c398cSBarry Smith . man - manual page with additional information on option 131294ae4db5SBarry Smith - currentvalue - the current value 131353acd3b1SBarry Smith 131453acd3b1SBarry Smith Output Parameter: 131553acd3b1SBarry Smith . flg - PETSC_TRUE or PETSC_FALSE 131653acd3b1SBarry Smith . set - PETSC_TRUE if found, else PETSC_FALSE 131753acd3b1SBarry Smith 1318*2efd9cb1SBarry Smith Notes: 1319*2efd9cb1SBarry Smith TRUE, true, YES, yes, nostring, and 1 all translate to PETSC_TRUE 1320*2efd9cb1SBarry Smith FALSE, false, NO, no, and 0 all translate to PETSC_FALSE 1321*2efd9cb1SBarry Smith 1322*2efd9cb1SBarry Smith If the option is given, but no value is provided, then flg and set are both given the value PETSC_TRUE. That is -requested_bool 1323*2efd9cb1SBarry Smith is equivalent to -requested_bool true 1324*2efd9cb1SBarry Smith 1325*2efd9cb1SBarry Smith If the user does not supply the option at all flg is NOT changed. Thus 1326*2efd9cb1SBarry Smith you should ALWAYS initialize the flg if you access it without first checking if the set flag is true. 1327*2efd9cb1SBarry Smith 132853acd3b1SBarry Smith Level: beginner 132953acd3b1SBarry Smith 133053acd3b1SBarry Smith Concepts: options database^logical 133153acd3b1SBarry Smith 133253acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 133353acd3b1SBarry Smith 1334c5929fdfSBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 1335acfcf0e5SJed Brown PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1336acfcf0e5SJed Brown PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 133753acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 133853acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1339acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1340a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 134153acd3b1SBarry Smith @*/ 13424416b707SBarry Smith PetscErrorCode PetscOptionsBool_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool *flg,PetscBool *set) 134353acd3b1SBarry Smith { 134453acd3b1SBarry Smith PetscErrorCode ierr; 1345ace3abfcSBarry Smith PetscBool iset; 13464416b707SBarry Smith PetscOptionItem amsopt; 134753acd3b1SBarry Smith 134853acd3b1SBarry Smith PetscFunctionBegin; 1349e55864a3SBarry Smith if (!PetscOptionsObject->count) { 13504416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr); 1351ace3abfcSBarry Smith ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr); 1352a297a907SKarl Rupp 135394ae4db5SBarry Smith *(PetscBool*)amsopt->data = currentvalue; 1354af6d86caSBarry Smith } 1355c5929fdfSBarry Smith ierr = PetscOptionsGetBool(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr); 135653acd3b1SBarry Smith if (set) *set = iset; 13571a1499c8SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 135894ae4db5SBarry Smith const char *v = PetscBools[currentvalue]; 13591a1499c8SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr); 136053acd3b1SBarry Smith } 136153acd3b1SBarry Smith PetscFunctionReturn(0); 136253acd3b1SBarry Smith } 136353acd3b1SBarry Smith 136453acd3b1SBarry Smith /*@C 136553acd3b1SBarry Smith PetscOptionsRealArray - Gets an array of double values for a particular 136653acd3b1SBarry Smith option in the database. The values must be separated with commas with 136753acd3b1SBarry Smith no intervening spaces. 136853acd3b1SBarry Smith 13693f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 137053acd3b1SBarry Smith 137153acd3b1SBarry Smith Input Parameters: 137253acd3b1SBarry Smith + opt - the option one is seeking 137353acd3b1SBarry Smith . text - short string describing option 137453acd3b1SBarry Smith . man - manual page for option 137553acd3b1SBarry Smith - nmax - maximum number of values 137653acd3b1SBarry Smith 137753acd3b1SBarry Smith Output Parameter: 137853acd3b1SBarry Smith + value - location to copy values 137953acd3b1SBarry Smith . nmax - actual number of values found 138053acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 138153acd3b1SBarry Smith 138253acd3b1SBarry Smith Level: beginner 138353acd3b1SBarry Smith 138453acd3b1SBarry Smith Notes: 138553acd3b1SBarry Smith The user should pass in an array of doubles 138653acd3b1SBarry Smith 138753acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 138853acd3b1SBarry Smith 138953acd3b1SBarry Smith Concepts: options database^array of strings 139053acd3b1SBarry Smith 1391c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1392acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 139353acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 139453acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1395acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1396a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 139753acd3b1SBarry Smith @*/ 13984416b707SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool *set) 139953acd3b1SBarry Smith { 140053acd3b1SBarry Smith PetscErrorCode ierr; 140153acd3b1SBarry Smith PetscInt i; 14024416b707SBarry Smith PetscOptionItem amsopt; 140353acd3b1SBarry Smith 140453acd3b1SBarry Smith PetscFunctionBegin; 1405e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1406e26ddf31SBarry Smith PetscReal *vals; 1407e26ddf31SBarry Smith 14084416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr); 1409e55864a3SBarry Smith ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr); 1410e26ddf31SBarry Smith vals = (PetscReal*)amsopt->data; 1411e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1412e26ddf31SBarry Smith amsopt->arraylength = *n; 1413e26ddf31SBarry Smith } 1414c5929fdfSBarry Smith ierr = PetscOptionsGetRealArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1415e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1416a519f713SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr); 141753acd3b1SBarry Smith for (i=1; i<*n; i++) { 1418e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr); 141953acd3b1SBarry Smith } 1420e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 142153acd3b1SBarry Smith } 142253acd3b1SBarry Smith PetscFunctionReturn(0); 142353acd3b1SBarry Smith } 142453acd3b1SBarry Smith 1425050cccc3SHong Zhang /*@C 1426050cccc3SHong Zhang PetscOptionsScalarArray - Gets an array of Scalar values for a particular 1427050cccc3SHong Zhang option in the database. The values must be separated with commas with 1428050cccc3SHong Zhang no intervening spaces. 1429050cccc3SHong Zhang 1430050cccc3SHong Zhang Logically Collective on the communicator passed in PetscOptionsBegin() 1431050cccc3SHong Zhang 1432050cccc3SHong Zhang Input Parameters: 1433050cccc3SHong Zhang + opt - the option one is seeking 1434050cccc3SHong Zhang . text - short string describing option 1435050cccc3SHong Zhang . man - manual page for option 1436050cccc3SHong Zhang - nmax - maximum number of values 1437050cccc3SHong Zhang 1438050cccc3SHong Zhang Output Parameter: 1439050cccc3SHong Zhang + value - location to copy values 1440050cccc3SHong Zhang . nmax - actual number of values found 1441050cccc3SHong Zhang - set - PETSC_TRUE if found, else PETSC_FALSE 1442050cccc3SHong Zhang 1443050cccc3SHong Zhang Level: beginner 1444050cccc3SHong Zhang 1445050cccc3SHong Zhang Notes: 1446050cccc3SHong Zhang The user should pass in an array of doubles 1447050cccc3SHong Zhang 1448050cccc3SHong Zhang Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1449050cccc3SHong Zhang 1450050cccc3SHong Zhang Concepts: options database^array of strings 1451050cccc3SHong Zhang 1452c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1453050cccc3SHong Zhang PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1454050cccc3SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1455050cccc3SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1456050cccc3SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1457050cccc3SHong Zhang PetscOptionsFList(), PetscOptionsEList() 1458050cccc3SHong Zhang @*/ 14594416b707SBarry Smith PetscErrorCode PetscOptionsScalarArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool *set) 1460050cccc3SHong Zhang { 1461050cccc3SHong Zhang PetscErrorCode ierr; 1462050cccc3SHong Zhang PetscInt i; 14634416b707SBarry Smith PetscOptionItem amsopt; 1464050cccc3SHong Zhang 1465050cccc3SHong Zhang PetscFunctionBegin; 1466050cccc3SHong Zhang if (!PetscOptionsObject->count) { 1467050cccc3SHong Zhang PetscScalar *vals; 1468050cccc3SHong Zhang 14694416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr); 1470050cccc3SHong Zhang ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr); 1471050cccc3SHong Zhang vals = (PetscScalar*)amsopt->data; 1472050cccc3SHong Zhang for (i=0; i<*n; i++) vals[i] = value[i]; 1473050cccc3SHong Zhang amsopt->arraylength = *n; 1474050cccc3SHong Zhang } 1475c5929fdfSBarry Smith ierr = PetscOptionsGetScalarArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1476050cccc3SHong Zhang if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 147795f3a755SJose 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); 1478050cccc3SHong Zhang for (i=1; i<*n; i++) { 147995f3a755SJose E. Roman ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g+%gi",(double)PetscRealPart(value[i]),(double)PetscImaginaryPart(value[i]));CHKERRQ(ierr); 1480050cccc3SHong Zhang } 1481050cccc3SHong Zhang ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 1482050cccc3SHong Zhang } 1483050cccc3SHong Zhang PetscFunctionReturn(0); 1484050cccc3SHong Zhang } 148553acd3b1SBarry Smith 148653acd3b1SBarry Smith /*@C 148753acd3b1SBarry Smith PetscOptionsIntArray - Gets an array of integers for a particular 1488b32a342fSShri Abhyankar option in the database. 148953acd3b1SBarry Smith 14903f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 149153acd3b1SBarry Smith 149253acd3b1SBarry Smith Input Parameters: 149353acd3b1SBarry Smith + opt - the option one is seeking 149453acd3b1SBarry Smith . text - short string describing option 149553acd3b1SBarry Smith . man - manual page for option 1496f8a50e2bSBarry Smith - n - maximum number of values 149753acd3b1SBarry Smith 149853acd3b1SBarry Smith Output Parameter: 149953acd3b1SBarry Smith + value - location to copy values 1500f8a50e2bSBarry Smith . n - actual number of values found 150153acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 150253acd3b1SBarry Smith 150353acd3b1SBarry Smith Level: beginner 150453acd3b1SBarry Smith 150553acd3b1SBarry Smith Notes: 1506b32a342fSShri Abhyankar The array can be passed as 1507bebe2cf6SSatish Balay a comma separated list: 0,1,2,3,4,5,6,7 15080fd488f5SShri Abhyankar a range (start-end+1): 0-8 15090fd488f5SShri Abhyankar a range with given increment (start-end+1:inc): 0-7:2 1510bebe2cf6SSatish Balay a combination of values and ranges separated by commas: 0,1-8,8-15:2 1511b32a342fSShri Abhyankar 1512b32a342fSShri Abhyankar There must be no intervening spaces between the values. 151353acd3b1SBarry Smith 151453acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 151553acd3b1SBarry Smith 1516b32a342fSShri Abhyankar Concepts: options database^array of ints 151753acd3b1SBarry Smith 1518c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1519acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 152053acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 152153acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1522acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1523a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray() 152453acd3b1SBarry Smith @*/ 15254416b707SBarry Smith PetscErrorCode PetscOptionsIntArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool *set) 152653acd3b1SBarry Smith { 152753acd3b1SBarry Smith PetscErrorCode ierr; 152853acd3b1SBarry Smith PetscInt i; 15294416b707SBarry Smith PetscOptionItem amsopt; 153053acd3b1SBarry Smith 153153acd3b1SBarry Smith PetscFunctionBegin; 1532e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1533e26ddf31SBarry Smith PetscInt *vals; 1534e26ddf31SBarry Smith 15354416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr); 1536854ce69bSBarry Smith ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr); 1537e26ddf31SBarry Smith vals = (PetscInt*)amsopt->data; 1538e26ddf31SBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 1539e26ddf31SBarry Smith amsopt->arraylength = *n; 1540e26ddf31SBarry Smith } 1541c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1542e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1543e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 154453acd3b1SBarry Smith for (i=1; i<*n; i++) { 1545e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 154653acd3b1SBarry Smith } 1547e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 154853acd3b1SBarry Smith } 154953acd3b1SBarry Smith PetscFunctionReturn(0); 155053acd3b1SBarry Smith } 155153acd3b1SBarry Smith 155253acd3b1SBarry Smith /*@C 155353acd3b1SBarry Smith PetscOptionsStringArray - Gets an array of string values for a particular 155453acd3b1SBarry Smith option in the database. The values must be separated with commas with 155553acd3b1SBarry Smith no intervening spaces. 155653acd3b1SBarry Smith 15573f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 155853acd3b1SBarry Smith 155953acd3b1SBarry Smith Input Parameters: 156053acd3b1SBarry Smith + opt - the option one is seeking 156153acd3b1SBarry Smith . text - short string describing option 156253acd3b1SBarry Smith . man - manual page for option 156353acd3b1SBarry Smith - nmax - maximum number of strings 156453acd3b1SBarry Smith 156553acd3b1SBarry Smith Output Parameter: 156653acd3b1SBarry Smith + value - location to copy strings 156753acd3b1SBarry Smith . nmax - actual number of strings found 156853acd3b1SBarry Smith - set - PETSC_TRUE if found, else PETSC_FALSE 156953acd3b1SBarry Smith 157053acd3b1SBarry Smith Level: beginner 157153acd3b1SBarry Smith 157253acd3b1SBarry Smith Notes: 157353acd3b1SBarry Smith The user should pass in an array of pointers to char, to hold all the 157453acd3b1SBarry Smith strings returned by this function. 157553acd3b1SBarry Smith 157653acd3b1SBarry Smith The user is responsible for deallocating the strings that are 157753acd3b1SBarry Smith returned. The Fortran interface for this routine is not supported. 157853acd3b1SBarry Smith 157953acd3b1SBarry Smith Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 158053acd3b1SBarry Smith 158153acd3b1SBarry Smith Concepts: options database^array of strings 158253acd3b1SBarry Smith 1583c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1584acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 158553acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 158653acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1587acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1588a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 158953acd3b1SBarry Smith @*/ 15904416b707SBarry Smith PetscErrorCode PetscOptionsStringArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool *set) 159153acd3b1SBarry Smith { 159253acd3b1SBarry Smith PetscErrorCode ierr; 15934416b707SBarry Smith PetscOptionItem amsopt; 159453acd3b1SBarry Smith 159553acd3b1SBarry Smith PetscFunctionBegin; 1596e55864a3SBarry Smith if (!PetscOptionsObject->count) { 15974416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr); 1598854ce69bSBarry Smith ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr); 1599a297a907SKarl Rupp 16001ae3d29cSBarry Smith amsopt->arraylength = *nmax; 16011ae3d29cSBarry Smith } 1602c5929fdfSBarry Smith ierr = PetscOptionsGetStringArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr); 1603e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1604e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr); 160553acd3b1SBarry Smith } 160653acd3b1SBarry Smith PetscFunctionReturn(0); 160753acd3b1SBarry Smith } 160853acd3b1SBarry Smith 1609e2446a98SMatthew Knepley /*@C 1610acfcf0e5SJed Brown PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular 1611e2446a98SMatthew Knepley option in the database. The values must be separated with commas with 1612e2446a98SMatthew Knepley no intervening spaces. 1613e2446a98SMatthew Knepley 16143f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 1615e2446a98SMatthew Knepley 1616e2446a98SMatthew Knepley Input Parameters: 1617e2446a98SMatthew Knepley + opt - the option one is seeking 1618e2446a98SMatthew Knepley . text - short string describing option 1619e2446a98SMatthew Knepley . man - manual page for option 1620e2446a98SMatthew Knepley - nmax - maximum number of values 1621e2446a98SMatthew Knepley 1622e2446a98SMatthew Knepley Output Parameter: 1623e2446a98SMatthew Knepley + value - location to copy values 1624e2446a98SMatthew Knepley . nmax - actual number of values found 1625e2446a98SMatthew Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 1626e2446a98SMatthew Knepley 1627e2446a98SMatthew Knepley Level: beginner 1628e2446a98SMatthew Knepley 1629e2446a98SMatthew Knepley Notes: 1630e2446a98SMatthew Knepley The user should pass in an array of doubles 1631e2446a98SMatthew Knepley 1632e2446a98SMatthew Knepley Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 1633e2446a98SMatthew Knepley 1634e2446a98SMatthew Knepley Concepts: options database^array of strings 1635e2446a98SMatthew Knepley 1636c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1637acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 1638e2446a98SMatthew Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1639e2446a98SMatthew Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1640acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1641a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 1642e2446a98SMatthew Knepley @*/ 16434416b707SBarry Smith PetscErrorCode PetscOptionsBoolArray_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set) 1644e2446a98SMatthew Knepley { 1645e2446a98SMatthew Knepley PetscErrorCode ierr; 1646e2446a98SMatthew Knepley PetscInt i; 16474416b707SBarry Smith PetscOptionItem amsopt; 1648e2446a98SMatthew Knepley 1649e2446a98SMatthew Knepley PetscFunctionBegin; 1650e55864a3SBarry Smith if (!PetscOptionsObject->count) { 1651ace3abfcSBarry Smith PetscBool *vals; 16521ae3d29cSBarry Smith 16534416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr); 16541a1499c8SBarry Smith ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr); 1655ace3abfcSBarry Smith vals = (PetscBool*)amsopt->data; 16561ae3d29cSBarry Smith for (i=0; i<*n; i++) vals[i] = value[i]; 16571ae3d29cSBarry Smith amsopt->arraylength = *n; 16581ae3d29cSBarry Smith } 1659c5929fdfSBarry Smith ierr = PetscOptionsGetBoolArray(PetscOptionsObject->options,PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr); 1660e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1661e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr); 1662e2446a98SMatthew Knepley for (i=1; i<*n; i++) { 1663e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr); 1664e2446a98SMatthew Knepley } 1665e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr); 1666e2446a98SMatthew Knepley } 1667e2446a98SMatthew Knepley PetscFunctionReturn(0); 1668e2446a98SMatthew Knepley } 1669e2446a98SMatthew Knepley 16708cc676e6SMatthew G Knepley /*@C 1671d1da0b69SBarry Smith PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user 16728cc676e6SMatthew G Knepley 16738cc676e6SMatthew G Knepley Logically Collective on the communicator passed in PetscOptionsBegin() 16748cc676e6SMatthew G Knepley 16758cc676e6SMatthew G Knepley Input Parameters: 16768cc676e6SMatthew G Knepley + opt - option name 16778cc676e6SMatthew G Knepley . text - short string that describes the option 16788cc676e6SMatthew G Knepley - man - manual page with additional information on option 16798cc676e6SMatthew G Knepley 16808cc676e6SMatthew G Knepley Output Parameter: 16818cc676e6SMatthew G Knepley + viewer - the viewer 16828cc676e6SMatthew G Knepley - set - PETSC_TRUE if found, else PETSC_FALSE 16838cc676e6SMatthew G Knepley 16848cc676e6SMatthew G Knepley Level: beginner 16858cc676e6SMatthew G Knepley 16868cc676e6SMatthew G Knepley Concepts: options database^has int 16878cc676e6SMatthew G Knepley 16888cc676e6SMatthew G Knepley Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 16898cc676e6SMatthew G Knepley 16905a7113b9SPatrick Sanan See PetscOptionsGetViewer() for the format of the supplied viewer and its options 16918cc676e6SMatthew G Knepley 1692c5929fdfSBarry Smith .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(NULL,), 16938cc676e6SMatthew G Knepley PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 16948cc676e6SMatthew G Knepley PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 16958cc676e6SMatthew G Knepley PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 16968cc676e6SMatthew G Knepley PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 16978cc676e6SMatthew G Knepley PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1698a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 16998cc676e6SMatthew G Knepley @*/ 17004416b707SBarry Smith PetscErrorCode PetscOptionsViewer_Private(PetscOptionItems *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool *set) 17018cc676e6SMatthew G Knepley { 17028cc676e6SMatthew G Knepley PetscErrorCode ierr; 17034416b707SBarry Smith PetscOptionItem amsopt; 17048cc676e6SMatthew G Knepley 17058cc676e6SMatthew G Knepley PetscFunctionBegin; 17061a1499c8SBarry Smith if (!PetscOptionsObject->count) { 17074416b707SBarry Smith ierr = PetscOptionItemCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr); 170864facd6cSBarry Smith /* must use system malloc since SAWs may free this */ 17095b02f95dSBarry Smith ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr); 17108cc676e6SMatthew G Knepley } 1711e55864a3SBarry Smith ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr); 1712e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1713e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr); 17148cc676e6SMatthew G Knepley } 17158cc676e6SMatthew G Knepley PetscFunctionReturn(0); 17168cc676e6SMatthew G Knepley } 17178cc676e6SMatthew G Knepley 171853acd3b1SBarry Smith 171953acd3b1SBarry Smith /*@C 1720b52f573bSBarry Smith PetscOptionsHead - Puts a heading before listing any more published options. Used, for example, 172153acd3b1SBarry Smith in KSPSetFromOptions_GMRES(). 172253acd3b1SBarry Smith 17233f9fe445SBarry Smith Logically Collective on the communicator passed in PetscOptionsBegin() 172453acd3b1SBarry Smith 172553acd3b1SBarry Smith Input Parameter: 172653acd3b1SBarry Smith . head - the heading text 172753acd3b1SBarry Smith 172853acd3b1SBarry Smith 172953acd3b1SBarry Smith Level: intermediate 173053acd3b1SBarry Smith 173153acd3b1SBarry Smith Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd() 173253acd3b1SBarry Smith 1733b52f573bSBarry Smith Can be followed by a call to PetscOptionsTail() in the same function. 173453acd3b1SBarry Smith 173553acd3b1SBarry Smith Concepts: options database^subheading 173653acd3b1SBarry Smith 1737c5929fdfSBarry Smith .seealso: PetscOptionsGetInt(NULL,), PetscOptionsGetReal(), 1738acfcf0e5SJed Brown PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(), 173953acd3b1SBarry Smith PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 174053acd3b1SBarry Smith PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1741acfcf0e5SJed Brown PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1742a264d7a6SBarry Smith PetscOptionsFList(), PetscOptionsEList() 174353acd3b1SBarry Smith @*/ 17444416b707SBarry Smith PetscErrorCode PetscOptionsHead(PetscOptionItems *PetscOptionsObject,const char head[]) 174553acd3b1SBarry Smith { 174653acd3b1SBarry Smith PetscErrorCode ierr; 174753acd3b1SBarry Smith 174853acd3b1SBarry Smith PetscFunctionBegin; 1749e55864a3SBarry Smith if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) { 1750e55864a3SBarry Smith ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s\n",head);CHKERRQ(ierr); 175153acd3b1SBarry Smith } 175253acd3b1SBarry Smith PetscFunctionReturn(0); 175353acd3b1SBarry Smith } 175453acd3b1SBarry Smith 175553acd3b1SBarry Smith 175653acd3b1SBarry Smith 175753acd3b1SBarry Smith 175853acd3b1SBarry Smith 175953acd3b1SBarry Smith 1760