xref: /petsc/src/sys/objects/aoptions.c (revision 88a9752c2b415d6528496bfe782c5e1be8388914)
17d0a6c19SBarry Smith 
253acd3b1SBarry Smith /*
33fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
43fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
553acd3b1SBarry Smith 
653acd3b1SBarry Smith */
753acd3b1SBarry Smith 
8afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
9665c2dedSJed Brown #include <petscviewer.h>
1053acd3b1SBarry Smith 
112aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
122aa6d131SJed Brown 
1353acd3b1SBarry Smith /*
1453acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
153fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1653acd3b1SBarry Smith 
1753acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1853acd3b1SBarry Smith */
19f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2053acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
2753acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
3253acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3353acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
346356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3753acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith 
410298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4253acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4361b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4553acd3b1SBarry Smith     }
4661b37b28SSatish Balay   }
4753acd3b1SBarry Smith   PetscFunctionReturn(0);
4853acd3b1SBarry Smith }
4953acd3b1SBarry Smith 
503194b578SJed Brown #undef __FUNCT__
513194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
553194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(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);
633194b578SJed Brown   PetscOptionsObject.object         = obj;
643194b578SJed Brown   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   }
723194b578SJed Brown   ierr = PetscOptionsBegin_Private(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 */
7953acd3b1SBarry Smith #undef __FUNCT__
8053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
81e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8253acd3b1SBarry Smith {
8353acd3b1SBarry Smith   int          ierr;
8453acd3b1SBarry Smith   PetscOptions next;
853be6e4c3SJed Brown   PetscBool    valid;
8653acd3b1SBarry Smith 
8753acd3b1SBarry Smith   PetscFunctionBegin;
883be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
893be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
903be6e4c3SJed Brown 
916e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9253acd3b1SBarry Smith   (*amsopt)->next = 0;
9353acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
946356e834SBarry Smith   (*amsopt)->type = t;
9553acd3b1SBarry Smith   (*amsopt)->data = 0;
9661b37b28SSatish Balay 
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9853acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
996356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10053acd3b1SBarry Smith 
101a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
102a297a907SKarl Rupp   else {
10353acd3b1SBarry Smith     next = PetscOptionsObject.next;
10453acd3b1SBarry Smith     while (next->next) next = next->next;
10553acd3b1SBarry Smith     next->next = *amsopt;
10653acd3b1SBarry Smith   }
10753acd3b1SBarry Smith   PetscFunctionReturn(0);
10853acd3b1SBarry Smith }
10953acd3b1SBarry Smith 
11053acd3b1SBarry Smith #undef __FUNCT__
111aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
112aee2cecaSBarry Smith /*
1133fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith     Collective on MPI_Comm
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith    Input Parameters:
1183fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1193fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1203fc1eb6aSBarry Smith -     str - location to store input
121aee2cecaSBarry Smith 
122aee2cecaSBarry Smith     Bugs:
123aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
124aee2cecaSBarry Smith 
125aee2cecaSBarry Smith */
1263fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
127aee2cecaSBarry Smith {
128330cf3c9SBarry Smith   size_t         i;
129aee2cecaSBarry Smith   char           c;
1303fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
131aee2cecaSBarry Smith   PetscErrorCode ierr;
132aee2cecaSBarry Smith 
133aee2cecaSBarry Smith   PetscFunctionBegin;
134aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
135aee2cecaSBarry Smith   if (!rank) {
136aee2cecaSBarry Smith     c = (char) getchar();
137aee2cecaSBarry Smith     i = 0;
138aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
139aee2cecaSBarry Smith       str[i++] = c;
140aee2cecaSBarry Smith       c = (char)getchar();
141aee2cecaSBarry Smith     }
142aee2cecaSBarry Smith     str[i] = 0;
143aee2cecaSBarry Smith   }
1444dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1453fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
146aee2cecaSBarry Smith   PetscFunctionReturn(0);
147aee2cecaSBarry Smith }
148aee2cecaSBarry Smith 
149aee2cecaSBarry Smith #undef __FUNCT__
150aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
151aee2cecaSBarry Smith /*
1523cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
153aee2cecaSBarry Smith 
154aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
155aee2cecaSBarry Smith 
1567781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1577781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1587781c08eSBarry Smith 
159aee2cecaSBarry Smith     Bugs:
1607781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1613cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
162aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
163aee2cecaSBarry Smith 
1643cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1653cc1e11dSBarry Smith      address space and communicating with the PETSc program
1663cc1e11dSBarry Smith 
167aee2cecaSBarry Smith */
168aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1696356e834SBarry Smith {
1706356e834SBarry Smith   PetscErrorCode ierr;
1716356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1726356e834SBarry Smith   char           str[512];
173a4404d99SBarry Smith   PetscInt       id;
1747781c08eSBarry Smith   PetscBool      bid;
175a4404d99SBarry Smith   PetscReal      ir,*valr;
176330cf3c9SBarry Smith   PetscInt       *vald;
177330cf3c9SBarry Smith   size_t         i;
1786356e834SBarry Smith 
179e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1806356e834SBarry Smith   while (next) {
1816356e834SBarry Smith     switch (next->type) {
1826356e834SBarry Smith     case OPTION_HEAD:
1836356e834SBarry Smith       break;
184e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
185e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
186e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
187e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
188e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
189e26ddf31SBarry Smith         if (i < next->arraylength-1) {
190e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
191e26ddf31SBarry Smith         }
192e26ddf31SBarry Smith       }
193e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
194e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
195e26ddf31SBarry Smith       if (str[0]) {
196e26ddf31SBarry Smith         PetscToken token;
197e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
198e26ddf31SBarry Smith         size_t     len;
199e26ddf31SBarry Smith         char       *value;
200ace3abfcSBarry Smith         PetscBool  foundrange;
201e26ddf31SBarry Smith 
202e26ddf31SBarry Smith         next->set = PETSC_TRUE;
203e26ddf31SBarry Smith         value     = str;
204e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
205e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
206e26ddf31SBarry Smith         while (n < nmax) {
207e26ddf31SBarry Smith           if (!value) break;
208e26ddf31SBarry Smith 
209e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
210e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
211e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
212e26ddf31SBarry Smith           if (value[0] == '-') i=2;
213e26ddf31SBarry Smith           else i=1;
214330cf3c9SBarry Smith           for (;i<len; i++) {
215e26ddf31SBarry Smith             if (value[i] == '-') {
216e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
217e26ddf31SBarry Smith               value[i] = 0;
218cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
219cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
220e32f2f54SBarry 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);
221e32f2f54SBarry 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);
222e26ddf31SBarry Smith               for (; start<end; start++) {
223e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
224e26ddf31SBarry Smith               }
225e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
226e26ddf31SBarry Smith               break;
227e26ddf31SBarry Smith             }
228e26ddf31SBarry Smith           }
229e26ddf31SBarry Smith           if (!foundrange) {
230cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
231e26ddf31SBarry Smith             dvalue++;
232e26ddf31SBarry Smith             n++;
233e26ddf31SBarry Smith           }
234e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
235e26ddf31SBarry Smith         }
2368c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
237e26ddf31SBarry Smith       }
238e26ddf31SBarry Smith       break;
239e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
240e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
241e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
242e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
243e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
244e26ddf31SBarry Smith         if (i < next->arraylength-1) {
245e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
246e26ddf31SBarry Smith         }
247e26ddf31SBarry Smith       }
248e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
249e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
250e26ddf31SBarry Smith       if (str[0]) {
251e26ddf31SBarry Smith         PetscToken token;
252e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
253e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
254e26ddf31SBarry Smith         char       *value;
255e26ddf31SBarry Smith 
256e26ddf31SBarry Smith         next->set = PETSC_TRUE;
257e26ddf31SBarry Smith         value     = str;
258e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
259e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
260e26ddf31SBarry Smith         while (n < nmax) {
261e26ddf31SBarry Smith           if (!value) break;
262cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
263e26ddf31SBarry Smith           dvalue++;
264e26ddf31SBarry Smith           n++;
265e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
266e26ddf31SBarry Smith         }
2678c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
268e26ddf31SBarry Smith       }
269e26ddf31SBarry Smith       break;
2706356e834SBarry Smith     case OPTION_INT:
271e26ddf31SBarry 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);
2723fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2733fc1eb6aSBarry Smith       if (str[0]) {
274c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
275c272547aSJed Brown         sscanf(str,"%lld",&id);
276c272547aSJed Brown #else
277aee2cecaSBarry Smith         sscanf(str,"%d",&id);
278c272547aSJed Brown #endif
279aee2cecaSBarry Smith         next->set = PETSC_TRUE;
280a297a907SKarl Rupp 
281aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
282aee2cecaSBarry Smith       }
283aee2cecaSBarry Smith       break;
284aee2cecaSBarry Smith     case OPTION_REAL:
285e26ddf31SBarry 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);
2863fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2873fc1eb6aSBarry Smith       if (str[0]) {
288ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
289a4404d99SBarry Smith         sscanf(str,"%e",&ir);
290ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
291aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
292ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
293d9822059SBarry Smith         ir = strtoflt128(str,0);
294d9822059SBarry Smith #else
295513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
296a4404d99SBarry Smith #endif
297aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
298aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
299aee2cecaSBarry Smith       }
300aee2cecaSBarry Smith       break;
3017781c08eSBarry Smith     case OPTION_BOOL:
3027781c08eSBarry 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);
3037781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3047781c08eSBarry Smith       if (str[0]) {
3057781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3067781c08eSBarry Smith         next->set = PETSC_TRUE;
3077781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3087781c08eSBarry Smith       }
3097781c08eSBarry Smith       break;
310aee2cecaSBarry Smith     case OPTION_STRING:
311e26ddf31SBarry 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);
3123fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3133fc1eb6aSBarry Smith       if (str[0]) {
314aee2cecaSBarry Smith         next->set = PETSC_TRUE;
3157781c08eSBarry Smith         next->data = (void*)strdup(str);
3166356e834SBarry Smith       }
3176356e834SBarry Smith       break;
3183cc1e11dSBarry Smith     case OPTION_LIST:
319140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3203cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3213cc1e11dSBarry Smith       if (str[0]) {
3223cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3233cc1e11dSBarry Smith         next->set = PETSC_TRUE;
3247781c08eSBarry Smith         next->data = (void*)strdup(str);
3253cc1e11dSBarry Smith       }
3263cc1e11dSBarry Smith       break;
327b432afdaSMatthew Knepley     default:
328b432afdaSMatthew Knepley       break;
3296356e834SBarry Smith     }
3306356e834SBarry Smith     next = next->next;
3316356e834SBarry Smith   }
3326356e834SBarry Smith   PetscFunctionReturn(0);
3336356e834SBarry Smith }
3346356e834SBarry Smith 
335e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
336e04113cfSBarry Smith #include <petscviewersaws.h>
337d5649816SBarry Smith 
338d5649816SBarry Smith static int count = 0;
339d5649816SBarry Smith 
340b3506946SBarry Smith #undef __FUNCT__
341e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
342e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
343d5649816SBarry Smith {
3442657e9d9SBarry Smith   PetscFunctionBegin;
345d5649816SBarry Smith   PetscFunctionReturn(0);
346d5649816SBarry Smith }
347d5649816SBarry Smith 
348d5649816SBarry Smith #undef __FUNCT__
3497781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
350b3506946SBarry Smith /*
3517781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
352b3506946SBarry Smith 
353b3506946SBarry Smith     Bugs:
354b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
355b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
356b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
357b3506946SBarry Smith 
358b3506946SBarry Smith 
359b3506946SBarry Smith */
360475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
361b3506946SBarry Smith {
362b3506946SBarry Smith   PetscErrorCode ierr;
363b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
364d5649816SBarry Smith   static int     mancount = 0;
365b3506946SBarry Smith   char           options[16];
3667aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
367*88a9752cSBarry Smith   char           manname[16],textname[16];
3682657e9d9SBarry Smith   char           dir[1024];
369b3506946SBarry Smith 
370b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
371b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
372a297a907SKarl Rupp 
373e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
3741bc75a8dSBarry Smith 
375d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
3767781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
3777781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
3782657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
3792657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
380b3506946SBarry Smith 
381b3506946SBarry Smith   while (next) {
382d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
3832657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
3847781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
385d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
3867781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
3877781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
3889f32e415SBarry Smith 
389b3506946SBarry Smith     switch (next->type) {
390b3506946SBarry Smith     case OPTION_HEAD:
391b3506946SBarry Smith       break;
392b3506946SBarry Smith     case OPTION_INT_ARRAY:
3937781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
3942657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
395b3506946SBarry Smith       break;
396b3506946SBarry Smith     case OPTION_REAL_ARRAY:
3977781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
3982657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
399b3506946SBarry Smith       break;
400b3506946SBarry Smith     case OPTION_INT:
4017781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4022657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
403b3506946SBarry Smith       break;
404b3506946SBarry Smith     case OPTION_REAL:
4057781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4062657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
407b3506946SBarry Smith       break;
4087781c08eSBarry Smith     case OPTION_BOOL:
4097781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4102657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4111ae3d29cSBarry Smith       break;
4127781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4137781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4142657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
41571f08665SBarry Smith       break;
416b3506946SBarry Smith     case OPTION_STRING:
4177781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4187781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4191ae3d29cSBarry Smith       break;
4201ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4217781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4222657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
423b3506946SBarry Smith       break;
424b3506946SBarry Smith     case OPTION_LIST:
4257781c08eSBarry Smith       {/*PetscInt ntext;*/
4267781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4277781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4287781c08eSBarry Smith       /*ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
4297781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4307781c08eSBarry Smith        PetscStackCallSAWs(SAWs_Register,(dir,next->edata,ntext-1,SAWs_READ,SAWs_STRING));*/
4311ae3d29cSBarry Smith       break;}
4321ae3d29cSBarry Smith     case OPTION_ELIST:
4337781c08eSBarry Smith       {/*PetscInt ntext = next->nlist; */
4347781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4357781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4367781c08eSBarry Smith       /*ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4371ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4382657e9d9SBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->text);CHKERRQ(ierr);
4397781c08eSBarry Smith        PetscStackCallSAWs(SAWs_Register,(dir,next->edata,ntext,SAWs_READ,SAWs_STRING));*/
44071f08665SBarry Smith       break;}
441b3506946SBarry Smith     default:
442b3506946SBarry Smith       break;
443b3506946SBarry Smith     }
444b3506946SBarry Smith     next = next->next;
445b3506946SBarry Smith   }
446b3506946SBarry Smith 
447b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4487aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
449b3506946SBarry Smith 
450*88a9752cSBarry Smith   /* determine if any values have been set in GUI */
451*88a9752cSBarry Smith   next = PetscOptionsObject.next;
452*88a9752cSBarry Smith   while (next) {
453*88a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
454*88a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
455*88a9752cSBarry Smith     printf("%d set \n",next->set);
456*88a9752cSBarry Smith     next = next->next;
457*88a9752cSBarry Smith   }
458*88a9752cSBarry Smith 
459b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
460b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
461b3506946SBarry Smith 
4629a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
463b3506946SBarry Smith   PetscFunctionReturn(0);
464b3506946SBarry Smith }
465b3506946SBarry Smith #endif
466b3506946SBarry Smith 
4676356e834SBarry Smith #undef __FUNCT__
46853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
46953acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47053acd3b1SBarry Smith {
47153acd3b1SBarry Smith   PetscErrorCode ierr;
4726356e834SBarry Smith   PetscOptions   last;
4736356e834SBarry Smith   char           option[256],value[1024],tmp[32];
474330cf3c9SBarry Smith   size_t         j;
47553acd3b1SBarry Smith 
47653acd3b1SBarry Smith   PetscFunctionBegin;
477aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
478b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
4797781c08eSBarry Smith #if defined(PETSC_HAVE_SAWS)
480475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
481b3506946SBarry Smith #else
48271f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
483b3506946SBarry Smith #endif
484aee2cecaSBarry Smith     }
485aee2cecaSBarry Smith   }
4866356e834SBarry Smith 
487c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
488c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4896356e834SBarry Smith 
4906356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4916356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49261b37b28SSatish Balay   /* reset alreadyprinted flag */
49361b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4943194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4950298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4966356e834SBarry Smith 
4976356e834SBarry Smith   while (PetscOptionsObject.next) {
4986356e834SBarry Smith     if (PetscOptionsObject.next->set) {
4996356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5006356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5016356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5026356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5036356e834SBarry Smith       } else {
5046356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5056356e834SBarry Smith       }
5066356e834SBarry Smith 
5076356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5086356e834SBarry Smith       case OPTION_HEAD:
5096356e834SBarry Smith         break;
510e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
511e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
512e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
513e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
514e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
515e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
516e26ddf31SBarry Smith         }
517e26ddf31SBarry Smith         break;
5186356e834SBarry Smith       case OPTION_INT:
5197a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5206356e834SBarry Smith         break;
5216356e834SBarry Smith       case OPTION_REAL:
5227a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5236356e834SBarry Smith         break;
5246356e834SBarry Smith       case OPTION_REAL_ARRAY:
5257a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5266356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5277a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5286356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5296356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5306356e834SBarry Smith         }
5316356e834SBarry Smith         break;
5327781c08eSBarry Smith       case OPTION_BOOL:
53371f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5346356e834SBarry Smith         break;
5357781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
536ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5371ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
538ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5391ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5401ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5411ae3d29cSBarry Smith         }
5421ae3d29cSBarry Smith         break;
5436356e834SBarry Smith       case OPTION_LIST:
5441ae3d29cSBarry Smith       case OPTION_ELIST:
5457781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5466356e834SBarry Smith         break;
5471ae3d29cSBarry Smith       case OPTION_STRING:
548475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5498da2146eSBarry Smith         break;
5501ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5511ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5521ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5531ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5541ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5551ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5561ae3d29cSBarry Smith         }
5576356e834SBarry Smith         break;
5586356e834SBarry Smith       }
5596356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5606356e834SBarry Smith     }
561503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
562503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5636356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56471f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
565a297a907SKarl Rupp 
5667781c08eSBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_LIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
5677781c08eSBarry Smith       free(PetscOptionsObject.next->data);
5687781c08eSBarry Smith     } else {
5697781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
5707781c08eSBarry Smith     }
5717781c08eSBarry Smith 
5726356e834SBarry Smith     last                    = PetscOptionsObject.next;
5736356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5746356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5756356e834SBarry Smith   }
5766356e834SBarry Smith   PetscOptionsObject.next = 0;
57753acd3b1SBarry Smith   PetscFunctionReturn(0);
57853acd3b1SBarry Smith }
57953acd3b1SBarry Smith 
58053acd3b1SBarry Smith #undef __FUNCT__
58153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
58253acd3b1SBarry Smith /*@C
58353acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58453acd3b1SBarry Smith 
5853f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58653acd3b1SBarry Smith 
58753acd3b1SBarry Smith    Input Parameters:
58853acd3b1SBarry Smith +  opt - option name
58953acd3b1SBarry Smith .  text - short string that describes the option
59053acd3b1SBarry Smith .  man - manual page with additional information on option
59153acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
59253acd3b1SBarry Smith -  defaultv - the default (current) value
59353acd3b1SBarry Smith 
59453acd3b1SBarry Smith    Output Parameter:
59553acd3b1SBarry Smith +  value - the  value to return
596b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
59753acd3b1SBarry Smith 
59853acd3b1SBarry Smith    Level: beginner
59953acd3b1SBarry Smith 
60053acd3b1SBarry Smith    Concepts: options database
60153acd3b1SBarry Smith 
60253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
60353acd3b1SBarry Smith 
60453acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60553acd3b1SBarry Smith 
60653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
607acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
608acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
61053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
611acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
61253acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
61353acd3b1SBarry Smith @*/
6147087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61553acd3b1SBarry Smith {
61653acd3b1SBarry Smith   PetscErrorCode ierr;
61753acd3b1SBarry Smith   PetscInt       ntext = 0;
618aa5bb8c0SSatish Balay   PetscInt       tval;
619ace3abfcSBarry Smith   PetscBool      tflg;
62053acd3b1SBarry Smith 
62153acd3b1SBarry Smith   PetscFunctionBegin;
62253acd3b1SBarry Smith   while (list[ntext++]) {
623e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62453acd3b1SBarry Smith   }
625e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62653acd3b1SBarry Smith   ntext -= 3;
627aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
628aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
629aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
630aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
63153acd3b1SBarry Smith   PetscFunctionReturn(0);
63253acd3b1SBarry Smith }
63353acd3b1SBarry Smith 
63453acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63553acd3b1SBarry Smith #undef __FUNCT__
63653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63753acd3b1SBarry Smith /*@C
63853acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63953acd3b1SBarry Smith 
6403f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64153acd3b1SBarry Smith 
64253acd3b1SBarry Smith    Input Parameters:
64353acd3b1SBarry Smith +  opt - option name
64453acd3b1SBarry Smith .  text - short string that describes the option
64553acd3b1SBarry Smith .  man - manual page with additional information on option
64653acd3b1SBarry Smith -  defaultv - the default (current) value
64753acd3b1SBarry Smith 
64853acd3b1SBarry Smith    Output Parameter:
64953acd3b1SBarry Smith +  value - the integer value to return
65053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith    Level: beginner
65353acd3b1SBarry Smith 
65453acd3b1SBarry Smith    Concepts: options database^has int
65553acd3b1SBarry Smith 
65653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65753acd3b1SBarry Smith 
65853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
659acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
660acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
663acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
66453acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
66553acd3b1SBarry Smith @*/
6667087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66753acd3b1SBarry Smith {
66853acd3b1SBarry Smith   PetscErrorCode ierr;
6696356e834SBarry Smith   PetscOptions   amsopt;
67053acd3b1SBarry Smith 
67153acd3b1SBarry Smith   PetscFunctionBegin;
672b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6736356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6746356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
675a297a907SKarl Rupp 
6766356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
677af6d86caSBarry Smith   }
67853acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6802aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
68153acd3b1SBarry Smith   }
68253acd3b1SBarry Smith   PetscFunctionReturn(0);
68353acd3b1SBarry Smith }
68453acd3b1SBarry Smith 
68553acd3b1SBarry Smith #undef __FUNCT__
68653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68753acd3b1SBarry Smith /*@C
68853acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68953acd3b1SBarry Smith 
6903f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith    Input Parameters:
69353acd3b1SBarry Smith +  opt - option name
69453acd3b1SBarry Smith .  text - short string that describes the option
69553acd3b1SBarry Smith .  man - manual page with additional information on option
696bcbf2dc5SJed Brown .  defaultv - the default (current) value
697bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69853acd3b1SBarry Smith 
69953acd3b1SBarry Smith    Output Parameter:
70053acd3b1SBarry Smith +  value - the value to return
70153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70253acd3b1SBarry Smith 
70353acd3b1SBarry Smith    Level: beginner
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith    Concepts: options database^has int
70653acd3b1SBarry Smith 
70753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70853acd3b1SBarry Smith 
7097fccdfe4SBarry 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).
7107fccdfe4SBarry Smith 
71153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
712acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
713acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
716acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71753acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71853acd3b1SBarry Smith @*/
7197087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
72053acd3b1SBarry Smith {
72153acd3b1SBarry Smith   PetscErrorCode ierr;
722aee2cecaSBarry Smith   PetscOptions   amsopt;
72353acd3b1SBarry Smith 
72453acd3b1SBarry Smith   PetscFunctionBegin;
725b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
726aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
7277781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
728af6d86caSBarry Smith   }
72953acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
73061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7312aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73253acd3b1SBarry Smith   }
73353acd3b1SBarry Smith   PetscFunctionReturn(0);
73453acd3b1SBarry Smith }
73553acd3b1SBarry Smith 
73653acd3b1SBarry Smith #undef __FUNCT__
73753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73853acd3b1SBarry Smith /*@C
73953acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
74053acd3b1SBarry Smith 
7413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74253acd3b1SBarry Smith 
74353acd3b1SBarry Smith    Input Parameters:
74453acd3b1SBarry Smith +  opt - option name
74553acd3b1SBarry Smith .  text - short string that describes the option
74653acd3b1SBarry Smith .  man - manual page with additional information on option
74753acd3b1SBarry Smith -  defaultv - the default (current) value
74853acd3b1SBarry Smith 
74953acd3b1SBarry Smith    Output Parameter:
75053acd3b1SBarry Smith +  value - the value to return
75153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75253acd3b1SBarry Smith 
75353acd3b1SBarry Smith    Level: beginner
75453acd3b1SBarry Smith 
75553acd3b1SBarry Smith    Concepts: options database^has int
75653acd3b1SBarry Smith 
75753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75853acd3b1SBarry Smith 
75953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
760acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
761acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
76253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
764acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
76553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
76653acd3b1SBarry Smith @*/
7677087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76853acd3b1SBarry Smith {
76953acd3b1SBarry Smith   PetscErrorCode ierr;
770538aa990SBarry Smith   PetscOptions   amsopt;
77153acd3b1SBarry Smith 
77253acd3b1SBarry Smith   PetscFunctionBegin;
773b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
774538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
775538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
776a297a907SKarl Rupp 
777538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
778538aa990SBarry Smith   }
77953acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
78061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7812aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78253acd3b1SBarry Smith   }
78353acd3b1SBarry Smith   PetscFunctionReturn(0);
78453acd3b1SBarry Smith }
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith #undef __FUNCT__
78753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78853acd3b1SBarry Smith /*@C
78953acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
79053acd3b1SBarry Smith 
7913f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79253acd3b1SBarry Smith 
79353acd3b1SBarry Smith    Input Parameters:
79453acd3b1SBarry Smith +  opt - option name
79553acd3b1SBarry Smith .  text - short string that describes the option
79653acd3b1SBarry Smith .  man - manual page with additional information on option
79753acd3b1SBarry Smith -  defaultv - the default (current) value
79853acd3b1SBarry Smith 
79953acd3b1SBarry Smith    Output Parameter:
80053acd3b1SBarry Smith +  value - the value to return
80153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith    Level: beginner
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith    Concepts: options database^has int
80653acd3b1SBarry Smith 
80753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80853acd3b1SBarry Smith 
80953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
810acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
811acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
814acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
81553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
81653acd3b1SBarry Smith @*/
8177087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81853acd3b1SBarry Smith {
81953acd3b1SBarry Smith   PetscErrorCode ierr;
82053acd3b1SBarry Smith 
82153acd3b1SBarry Smith   PetscFunctionBegin;
82253acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
82353acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82453acd3b1SBarry Smith #else
82553acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82653acd3b1SBarry Smith #endif
82753acd3b1SBarry Smith   PetscFunctionReturn(0);
82853acd3b1SBarry Smith }
82953acd3b1SBarry Smith 
83053acd3b1SBarry Smith #undef __FUNCT__
83153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
83253acd3b1SBarry Smith /*@C
83390d69ab7SBarry 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
83490d69ab7SBarry Smith                       its value is set to false.
83553acd3b1SBarry Smith 
8363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83753acd3b1SBarry Smith 
83853acd3b1SBarry Smith    Input Parameters:
83953acd3b1SBarry Smith +  opt - option name
84053acd3b1SBarry Smith .  text - short string that describes the option
84153acd3b1SBarry Smith -  man - manual page with additional information on option
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith    Output Parameter:
84453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith    Level: beginner
84753acd3b1SBarry Smith 
84853acd3b1SBarry Smith    Concepts: options database^has int
84953acd3b1SBarry Smith 
85053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85153acd3b1SBarry Smith 
85253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
853acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
854acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
857acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
85853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
85953acd3b1SBarry Smith @*/
8607087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
86153acd3b1SBarry Smith {
86253acd3b1SBarry Smith   PetscErrorCode ierr;
8631ae3d29cSBarry Smith   PetscOptions   amsopt;
86453acd3b1SBarry Smith 
86553acd3b1SBarry Smith   PetscFunctionBegin;
8661ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8677781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
868ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
869a297a907SKarl Rupp 
870ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8711ae3d29cSBarry Smith   }
87253acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
87361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8742aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
87553acd3b1SBarry Smith   }
87653acd3b1SBarry Smith   PetscFunctionReturn(0);
87753acd3b1SBarry Smith }
87853acd3b1SBarry Smith 
87953acd3b1SBarry Smith #undef __FUNCT__
88053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
88153acd3b1SBarry Smith /*@C
88253acd3b1SBarry Smith      PetscOptionsList - Puts a list of option values that a single one may be selected from
88353acd3b1SBarry Smith 
8843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88553acd3b1SBarry Smith 
88653acd3b1SBarry Smith    Input Parameters:
88753acd3b1SBarry Smith +  opt - option name
88853acd3b1SBarry Smith .  text - short string that describes the option
88953acd3b1SBarry Smith .  man - manual page with additional information on option
89053acd3b1SBarry Smith .  list - the possible choices
8913cc1e11dSBarry Smith .  defaultv - the default (current) value
8923cc1e11dSBarry Smith -  len - the length of the character array value
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith    Output Parameter:
89553acd3b1SBarry Smith +  value - the value to return
89653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    Level: intermediate
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
90353acd3b1SBarry Smith 
90453acd3b1SBarry Smith    To get a listing of all currently specified options,
90588c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
90653acd3b1SBarry Smith 
907eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
908eabe10d7SBarry Smith 
90953acd3b1SBarry Smith    Concepts: options database^list
91053acd3b1SBarry Smith 
91153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
912acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
91353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
915acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
916171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
91753acd3b1SBarry Smith @*/
918140e18c1SBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91953acd3b1SBarry Smith {
92053acd3b1SBarry Smith   PetscErrorCode ierr;
9213cc1e11dSBarry Smith   PetscOptions   amsopt;
92253acd3b1SBarry Smith 
92353acd3b1SBarry Smith   PetscFunctionBegin;
924b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9253cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
9267781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9273cc1e11dSBarry Smith     amsopt->flist = list;
9283cc1e11dSBarry Smith   }
92953acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
93061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
931140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
93253acd3b1SBarry Smith   }
93353acd3b1SBarry Smith   PetscFunctionReturn(0);
93453acd3b1SBarry Smith }
93553acd3b1SBarry Smith 
93653acd3b1SBarry Smith #undef __FUNCT__
93753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
93853acd3b1SBarry Smith /*@C
93953acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
94053acd3b1SBarry Smith 
9413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94253acd3b1SBarry Smith 
94353acd3b1SBarry Smith    Input Parameters:
94453acd3b1SBarry Smith +  opt - option name
94553acd3b1SBarry Smith .  ltext - short string that describes the option
94653acd3b1SBarry Smith .  man - manual page with additional information on option
94753acd3b1SBarry Smith .  list - the possible choices
94853acd3b1SBarry Smith .  ntext - number of choices
94953acd3b1SBarry Smith -  defaultv - the default (current) value
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith    Output Parameter:
95253acd3b1SBarry Smith +  value - the index of the value to return
95353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95453acd3b1SBarry Smith 
95553acd3b1SBarry Smith    Level: intermediate
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95853acd3b1SBarry Smith 
959140e18c1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFunctionList()
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith    Concepts: options database^list
96253acd3b1SBarry Smith 
96353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
964acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
967acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
968171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
96953acd3b1SBarry Smith @*/
9707087cfbeSBarry Smith PetscErrorCode  PetscOptionsEList(const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char defaultv[],PetscInt *value,PetscBool  *set)
97153acd3b1SBarry Smith {
97253acd3b1SBarry Smith   PetscErrorCode ierr;
97353acd3b1SBarry Smith   PetscInt       i;
9741ae3d29cSBarry Smith   PetscOptions   amsopt;
97553acd3b1SBarry Smith 
97653acd3b1SBarry Smith   PetscFunctionBegin;
9771ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
978d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9797781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9801ae3d29cSBarry Smith     amsopt->list  = list;
9811ae3d29cSBarry Smith     amsopt->nlist = ntext;
9821ae3d29cSBarry Smith   }
98353acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
98553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
98653acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
98753acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
98853acd3b1SBarry Smith     }
9892aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
99053acd3b1SBarry Smith   }
99153acd3b1SBarry Smith   PetscFunctionReturn(0);
99253acd3b1SBarry Smith }
99353acd3b1SBarry Smith 
99453acd3b1SBarry Smith #undef __FUNCT__
995acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
99653acd3b1SBarry Smith /*@C
997acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
998d5649816SBarry Smith        which at most a single value can be true.
99953acd3b1SBarry Smith 
10003f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100153acd3b1SBarry Smith 
100253acd3b1SBarry Smith    Input Parameters:
100353acd3b1SBarry Smith +  opt - option name
100453acd3b1SBarry Smith .  text - short string that describes the option
100553acd3b1SBarry Smith -  man - manual page with additional information on option
100653acd3b1SBarry Smith 
100753acd3b1SBarry Smith    Output Parameter:
100853acd3b1SBarry Smith .  flg - whether that option was set or not
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith    Level: intermediate
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101353acd3b1SBarry Smith 
1014acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
101553acd3b1SBarry Smith 
101653acd3b1SBarry Smith     Concepts: options database^logical group
101753acd3b1SBarry Smith 
101853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1019acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1022acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
102353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
102453acd3b1SBarry Smith @*/
10257087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
102653acd3b1SBarry Smith {
102753acd3b1SBarry Smith   PetscErrorCode ierr;
10281ae3d29cSBarry Smith   PetscOptions   amsopt;
102953acd3b1SBarry Smith 
103053acd3b1SBarry Smith   PetscFunctionBegin;
10311ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10327781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1033ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1034a297a907SKarl Rupp 
1035ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10361ae3d29cSBarry Smith   }
103768b16fdaSBarry Smith   *flg = PETSC_FALSE;
10380298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
103961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104053acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10412aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104253acd3b1SBarry Smith   }
104353acd3b1SBarry Smith   PetscFunctionReturn(0);
104453acd3b1SBarry Smith }
104553acd3b1SBarry Smith 
104653acd3b1SBarry Smith #undef __FUNCT__
1047acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
104853acd3b1SBarry Smith /*@C
1049acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1050d5649816SBarry Smith        which at most a single value can be true.
105153acd3b1SBarry Smith 
10523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105353acd3b1SBarry Smith 
105453acd3b1SBarry Smith    Input Parameters:
105553acd3b1SBarry Smith +  opt - option name
105653acd3b1SBarry Smith .  text - short string that describes the option
105753acd3b1SBarry Smith -  man - manual page with additional information on option
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith    Output Parameter:
106053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106153acd3b1SBarry Smith 
106253acd3b1SBarry Smith    Level: intermediate
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106553acd3b1SBarry Smith 
1066acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
106753acd3b1SBarry Smith 
106853acd3b1SBarry Smith     Concepts: options database^logical group
106953acd3b1SBarry Smith 
107053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1071acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1074acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
107553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
107653acd3b1SBarry Smith @*/
10777087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
107853acd3b1SBarry Smith {
107953acd3b1SBarry Smith   PetscErrorCode ierr;
10801ae3d29cSBarry Smith   PetscOptions   amsopt;
108153acd3b1SBarry Smith 
108253acd3b1SBarry Smith   PetscFunctionBegin;
10831ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10847781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1085ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1086a297a907SKarl Rupp 
1087ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10881ae3d29cSBarry Smith   }
108917326d04SJed Brown   *flg = PETSC_FALSE;
10900298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10922aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109353acd3b1SBarry Smith   }
109453acd3b1SBarry Smith   PetscFunctionReturn(0);
109553acd3b1SBarry Smith }
109653acd3b1SBarry Smith 
109753acd3b1SBarry Smith #undef __FUNCT__
1098acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
109953acd3b1SBarry Smith /*@C
1100acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1101d5649816SBarry Smith        which at most a single value can be true.
110253acd3b1SBarry Smith 
11033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110453acd3b1SBarry Smith 
110553acd3b1SBarry Smith    Input Parameters:
110653acd3b1SBarry Smith +  opt - option name
110753acd3b1SBarry Smith .  text - short string that describes the option
110853acd3b1SBarry Smith -  man - manual page with additional information on option
110953acd3b1SBarry Smith 
111053acd3b1SBarry Smith    Output Parameter:
111153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Level: intermediate
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111653acd3b1SBarry Smith 
1117acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith     Concepts: options database^logical group
112053acd3b1SBarry Smith 
112153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1122acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1125acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
112653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
112753acd3b1SBarry Smith @*/
11287087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
112953acd3b1SBarry Smith {
113053acd3b1SBarry Smith   PetscErrorCode ierr;
11311ae3d29cSBarry Smith   PetscOptions   amsopt;
113253acd3b1SBarry Smith 
113353acd3b1SBarry Smith   PetscFunctionBegin;
11341ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11357781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1136ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1137a297a907SKarl Rupp 
1138ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11391ae3d29cSBarry Smith   }
114017326d04SJed Brown   *flg = PETSC_FALSE;
11410298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11432aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114453acd3b1SBarry Smith   }
114553acd3b1SBarry Smith   PetscFunctionReturn(0);
114653acd3b1SBarry Smith }
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith #undef __FUNCT__
1149acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
115053acd3b1SBarry Smith /*@C
1151acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115253acd3b1SBarry Smith 
11533f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115453acd3b1SBarry Smith 
115553acd3b1SBarry Smith    Input Parameters:
115653acd3b1SBarry Smith +  opt - option name
115753acd3b1SBarry Smith .  text - short string that describes the option
115853acd3b1SBarry Smith -  man - manual page with additional information on option
115953acd3b1SBarry Smith 
116053acd3b1SBarry Smith    Output Parameter:
116153acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116253acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith    Level: beginner
116553acd3b1SBarry Smith 
116653acd3b1SBarry Smith    Concepts: options database^logical
116753acd3b1SBarry Smith 
116853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116953acd3b1SBarry Smith 
117053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1171acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1172acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
117353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1175acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
117653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
117753acd3b1SBarry Smith @*/
11787087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
117953acd3b1SBarry Smith {
118053acd3b1SBarry Smith   PetscErrorCode ierr;
1181ace3abfcSBarry Smith   PetscBool      iset;
1182aee2cecaSBarry Smith   PetscOptions   amsopt;
118353acd3b1SBarry Smith 
118453acd3b1SBarry Smith   PetscFunctionBegin;
1185b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
11867781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1187ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1188a297a907SKarl Rupp 
1189ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1190af6d86caSBarry Smith   }
1191acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119253acd3b1SBarry Smith   if (!iset) {
119353acd3b1SBarry Smith     if (flg) *flg = deflt;
119453acd3b1SBarry Smith   }
119553acd3b1SBarry Smith   if (set) *set = iset;
119661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1197ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
11982aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
119953acd3b1SBarry Smith   }
120053acd3b1SBarry Smith   PetscFunctionReturn(0);
120153acd3b1SBarry Smith }
120253acd3b1SBarry Smith 
120353acd3b1SBarry Smith #undef __FUNCT__
120453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
120553acd3b1SBarry Smith /*@C
120653acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
120753acd3b1SBarry Smith    option in the database. The values must be separated with commas with
120853acd3b1SBarry Smith    no intervening spaces.
120953acd3b1SBarry Smith 
12103f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121153acd3b1SBarry Smith 
121253acd3b1SBarry Smith    Input Parameters:
121353acd3b1SBarry Smith +  opt - the option one is seeking
121453acd3b1SBarry Smith .  text - short string describing option
121553acd3b1SBarry Smith .  man - manual page for option
121653acd3b1SBarry Smith -  nmax - maximum number of values
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith    Output Parameter:
121953acd3b1SBarry Smith +  value - location to copy values
122053acd3b1SBarry Smith .  nmax - actual number of values found
122153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Level: beginner
122453acd3b1SBarry Smith 
122553acd3b1SBarry Smith    Notes:
122653acd3b1SBarry Smith    The user should pass in an array of doubles
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122953acd3b1SBarry Smith 
123053acd3b1SBarry Smith    Concepts: options database^array of strings
123153acd3b1SBarry Smith 
123253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1233acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1236acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
123753acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
123853acd3b1SBarry Smith @*/
12397087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
124053acd3b1SBarry Smith {
124153acd3b1SBarry Smith   PetscErrorCode ierr;
124253acd3b1SBarry Smith   PetscInt       i;
1243e26ddf31SBarry Smith   PetscOptions   amsopt;
124453acd3b1SBarry Smith 
124553acd3b1SBarry Smith   PetscFunctionBegin;
1246b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1247e26ddf31SBarry Smith     PetscReal *vals;
1248e26ddf31SBarry Smith 
1249e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1250e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1251e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1252e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1253e26ddf31SBarry Smith     amsopt->arraylength = *n;
1254e26ddf31SBarry Smith   }
125553acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
125661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1257a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
125853acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1259a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
126053acd3b1SBarry Smith     }
12612aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126253acd3b1SBarry Smith   }
126353acd3b1SBarry Smith   PetscFunctionReturn(0);
126453acd3b1SBarry Smith }
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith 
126753acd3b1SBarry Smith #undef __FUNCT__
126853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
126953acd3b1SBarry Smith /*@C
127053acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1271b32a342fSShri Abhyankar    option in the database.
127253acd3b1SBarry Smith 
12733f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127453acd3b1SBarry Smith 
127553acd3b1SBarry Smith    Input Parameters:
127653acd3b1SBarry Smith +  opt - the option one is seeking
127753acd3b1SBarry Smith .  text - short string describing option
127853acd3b1SBarry Smith .  man - manual page for option
1279f8a50e2bSBarry Smith -  n - maximum number of values
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith    Output Parameter:
128253acd3b1SBarry Smith +  value - location to copy values
1283f8a50e2bSBarry Smith .  n - actual number of values found
128453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128553acd3b1SBarry Smith 
128653acd3b1SBarry Smith    Level: beginner
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith    Notes:
1289b32a342fSShri Abhyankar    The array can be passed as
1290b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12910fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12920fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12930fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1294b32a342fSShri Abhyankar 
1295b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
129653acd3b1SBarry Smith 
129753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129853acd3b1SBarry Smith 
1299b32a342fSShri Abhyankar    Concepts: options database^array of ints
130053acd3b1SBarry Smith 
130153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1302acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
130353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1305acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
130653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
130753acd3b1SBarry Smith @*/
13087087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
130953acd3b1SBarry Smith {
131053acd3b1SBarry Smith   PetscErrorCode ierr;
131153acd3b1SBarry Smith   PetscInt       i;
1312e26ddf31SBarry Smith   PetscOptions   amsopt;
131353acd3b1SBarry Smith 
131453acd3b1SBarry Smith   PetscFunctionBegin;
1315b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1316e26ddf31SBarry Smith     PetscInt *vals;
1317e26ddf31SBarry Smith 
1318e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1319e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1320e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1321e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1322e26ddf31SBarry Smith     amsopt->arraylength = *n;
1323e26ddf31SBarry Smith   }
132453acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
132561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132653acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
132753acd3b1SBarry Smith     for (i=1; i<*n; i++) {
132853acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
132953acd3b1SBarry Smith     }
13302aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133153acd3b1SBarry Smith   }
133253acd3b1SBarry Smith   PetscFunctionReturn(0);
133353acd3b1SBarry Smith }
133453acd3b1SBarry Smith 
133553acd3b1SBarry Smith #undef __FUNCT__
133653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
133753acd3b1SBarry Smith /*@C
133853acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
133953acd3b1SBarry Smith    option in the database. The values must be separated with commas with
134053acd3b1SBarry Smith    no intervening spaces.
134153acd3b1SBarry Smith 
13423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
134353acd3b1SBarry Smith 
134453acd3b1SBarry Smith    Input Parameters:
134553acd3b1SBarry Smith +  opt - the option one is seeking
134653acd3b1SBarry Smith .  text - short string describing option
134753acd3b1SBarry Smith .  man - manual page for option
134853acd3b1SBarry Smith -  nmax - maximum number of strings
134953acd3b1SBarry Smith 
135053acd3b1SBarry Smith    Output Parameter:
135153acd3b1SBarry Smith +  value - location to copy strings
135253acd3b1SBarry Smith .  nmax - actual number of strings found
135353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135453acd3b1SBarry Smith 
135553acd3b1SBarry Smith    Level: beginner
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith    Notes:
135853acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
135953acd3b1SBarry Smith    strings returned by this function.
136053acd3b1SBarry Smith 
136153acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136253acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
136353acd3b1SBarry Smith 
136453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136553acd3b1SBarry Smith 
136653acd3b1SBarry Smith    Concepts: options database^array of strings
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1369acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1372acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
137353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
137453acd3b1SBarry Smith @*/
13757087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
137653acd3b1SBarry Smith {
137753acd3b1SBarry Smith   PetscErrorCode ierr;
13781ae3d29cSBarry Smith   PetscOptions   amsopt;
137953acd3b1SBarry Smith 
138053acd3b1SBarry Smith   PetscFunctionBegin;
13811ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13821ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13831ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1384a297a907SKarl Rupp 
13851ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13861ae3d29cSBarry Smith   }
138753acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
138861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13892aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139053acd3b1SBarry Smith   }
139153acd3b1SBarry Smith   PetscFunctionReturn(0);
139253acd3b1SBarry Smith }
139353acd3b1SBarry Smith 
1394e2446a98SMatthew Knepley #undef __FUNCT__
1395acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1396e2446a98SMatthew Knepley /*@C
1397acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1398e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1399e2446a98SMatthew Knepley    no intervening spaces.
1400e2446a98SMatthew Knepley 
14013f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1402e2446a98SMatthew Knepley 
1403e2446a98SMatthew Knepley    Input Parameters:
1404e2446a98SMatthew Knepley +  opt - the option one is seeking
1405e2446a98SMatthew Knepley .  text - short string describing option
1406e2446a98SMatthew Knepley .  man - manual page for option
1407e2446a98SMatthew Knepley -  nmax - maximum number of values
1408e2446a98SMatthew Knepley 
1409e2446a98SMatthew Knepley    Output Parameter:
1410e2446a98SMatthew Knepley +  value - location to copy values
1411e2446a98SMatthew Knepley .  nmax - actual number of values found
1412e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1413e2446a98SMatthew Knepley 
1414e2446a98SMatthew Knepley    Level: beginner
1415e2446a98SMatthew Knepley 
1416e2446a98SMatthew Knepley    Notes:
1417e2446a98SMatthew Knepley    The user should pass in an array of doubles
1418e2446a98SMatthew Knepley 
1419e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1420e2446a98SMatthew Knepley 
1421e2446a98SMatthew Knepley    Concepts: options database^array of strings
1422e2446a98SMatthew Knepley 
1423e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1424acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1425e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1426e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1427acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1428e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1429e2446a98SMatthew Knepley @*/
14307087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1431e2446a98SMatthew Knepley {
1432e2446a98SMatthew Knepley   PetscErrorCode ierr;
1433e2446a98SMatthew Knepley   PetscInt       i;
14341ae3d29cSBarry Smith   PetscOptions   amsopt;
1435e2446a98SMatthew Knepley 
1436e2446a98SMatthew Knepley   PetscFunctionBegin;
14371ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1438ace3abfcSBarry Smith     PetscBool *vals;
14391ae3d29cSBarry Smith 
14407781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1441ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1442ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14431ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14441ae3d29cSBarry Smith     amsopt->arraylength = *n;
14451ae3d29cSBarry Smith   }
1446acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1447e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1448e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1449e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1450e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1451e2446a98SMatthew Knepley     }
14522aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1453e2446a98SMatthew Knepley   }
1454e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1455e2446a98SMatthew Knepley }
1456e2446a98SMatthew Knepley 
14578cc676e6SMatthew G Knepley #undef __FUNCT__
14588cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14598cc676e6SMatthew G Knepley /*@C
14608cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14618cc676e6SMatthew G Knepley 
14628cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14638cc676e6SMatthew G Knepley 
14648cc676e6SMatthew G Knepley    Input Parameters:
14658cc676e6SMatthew G Knepley +  opt - option name
14668cc676e6SMatthew G Knepley .  text - short string that describes the option
14678cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14688cc676e6SMatthew G Knepley 
14698cc676e6SMatthew G Knepley    Output Parameter:
14708cc676e6SMatthew G Knepley +  viewer - the viewer
14718cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14728cc676e6SMatthew G Knepley 
14738cc676e6SMatthew G Knepley    Level: beginner
14748cc676e6SMatthew G Knepley 
14758cc676e6SMatthew G Knepley    Concepts: options database^has int
14768cc676e6SMatthew G Knepley 
14778cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14788cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14798cc676e6SMatthew G Knepley $       ascii[:[filename][:format]]   defaults to stdout - format can be one of info, info_detailed, or matlab, for example ascii::info prints just the info
14808cc676e6SMatthew G Knepley $                                     about the object to standard out
14818cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14828cc676e6SMatthew G Knepley $       draw
14838cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14848cc676e6SMatthew G Knepley 
1485cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14868cc676e6SMatthew G Knepley 
14878cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14888cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14898cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14908cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14918cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14928cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14938cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14948cc676e6SMatthew G Knepley @*/
1495cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14968cc676e6SMatthew G Knepley {
14978cc676e6SMatthew G Knepley   PetscErrorCode ierr;
14988cc676e6SMatthew G Knepley   PetscOptions   amsopt;
14998cc676e6SMatthew G Knepley 
15008cc676e6SMatthew G Knepley   PetscFunctionBegin;
15018cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15028cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
15037781c08eSBarry Smith     amsopt->data = (void*)strdup("");
15048cc676e6SMatthew G Knepley   }
1505cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15068cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15078cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15088cc676e6SMatthew G Knepley   }
15098cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15108cc676e6SMatthew G Knepley }
15118cc676e6SMatthew G Knepley 
151253acd3b1SBarry Smith 
151353acd3b1SBarry Smith #undef __FUNCT__
151453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
151553acd3b1SBarry Smith /*@C
1516b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
151753acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
151853acd3b1SBarry Smith 
15193f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith    Input Parameter:
152253acd3b1SBarry Smith .   head - the heading text
152353acd3b1SBarry Smith 
152453acd3b1SBarry Smith 
152553acd3b1SBarry Smith    Level: intermediate
152653acd3b1SBarry Smith 
152753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152853acd3b1SBarry Smith 
1529b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
153053acd3b1SBarry Smith 
153153acd3b1SBarry Smith    Concepts: options database^subheading
153253acd3b1SBarry Smith 
153353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1534acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1537acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
153853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
153953acd3b1SBarry Smith @*/
15407087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154153acd3b1SBarry Smith {
154253acd3b1SBarry Smith   PetscErrorCode ierr;
154353acd3b1SBarry Smith 
154453acd3b1SBarry Smith   PetscFunctionBegin;
154561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
154653acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
154753acd3b1SBarry Smith   }
154853acd3b1SBarry Smith   PetscFunctionReturn(0);
154953acd3b1SBarry Smith }
155053acd3b1SBarry Smith 
155153acd3b1SBarry Smith 
155253acd3b1SBarry Smith 
155353acd3b1SBarry Smith 
155453acd3b1SBarry Smith 
155553acd3b1SBarry Smith 
1556