xref: /petsc/src/sys/objects/aoptions.c (revision 77b82169ce741be1f9ad766c5deaa64f2b843e90)
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;
315*77b82169SBarry Smith         ierr = PetscStrallocpy(str,(char**) &next->data);CHKERRQ(ierr);
3166356e834SBarry Smith       }
3176356e834SBarry Smith       break;
318a264d7a6SBarry Smith     case OPTION_FLIST:
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;
324*77b82169SBarry Smith         ierr = PetscStrallocpy(str,(char**) &next->data);CHKERRQ(ierr);
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;
36788a9752cSBarry 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;
424a264d7a6SBarry Smith     case OPTION_FLIST:
425a264d7a6SBarry Smith       {
426a264d7a6SBarry Smith       PetscInt ntext;
4277781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4287781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
429a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
430a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
431a264d7a6SBarry Smith       }
432a264d7a6SBarry Smith       break;
4331ae3d29cSBarry Smith     case OPTION_ELIST:
434a264d7a6SBarry Smith       {
435a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4367781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4377781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
438a264d7a6SBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char*),&next->edata);CHKERRQ(ierr);
4391ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
440a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
441a264d7a6SBarry Smith       }
442a264d7a6SBarry Smith       break;
443b3506946SBarry Smith     default:
444b3506946SBarry Smith       break;
445b3506946SBarry Smith     }
446b3506946SBarry Smith     next = next->next;
447b3506946SBarry Smith   }
448b3506946SBarry Smith 
449b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4507aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
451b3506946SBarry Smith 
45288a9752cSBarry Smith   /* determine if any values have been set in GUI */
45388a9752cSBarry Smith   next = PetscOptionsObject.next;
45488a9752cSBarry Smith   while (next) {
45588a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
45688a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
45788a9752cSBarry Smith     next = next->next;
45888a9752cSBarry Smith   }
45988a9752cSBarry Smith 
460b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
461b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
462b3506946SBarry Smith 
4639a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
464b3506946SBarry Smith   PetscFunctionReturn(0);
465b3506946SBarry Smith }
466b3506946SBarry Smith #endif
467b3506946SBarry Smith 
4686356e834SBarry Smith #undef __FUNCT__
46953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
47053acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47153acd3b1SBarry Smith {
47253acd3b1SBarry Smith   PetscErrorCode ierr;
4736356e834SBarry Smith   PetscOptions   last;
4746356e834SBarry Smith   char           option[256],value[1024],tmp[32];
475330cf3c9SBarry Smith   size_t         j;
47653acd3b1SBarry Smith 
47753acd3b1SBarry Smith   PetscFunctionBegin;
478aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
479b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
480a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
481475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
482b3506946SBarry Smith #else
48371f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
484b3506946SBarry Smith #endif
485aee2cecaSBarry Smith     }
486aee2cecaSBarry Smith   }
4876356e834SBarry Smith 
488c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
489c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4906356e834SBarry Smith 
4916356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4926356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49361b37b28SSatish Balay   /* reset alreadyprinted flag */
49461b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4953194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4960298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4976356e834SBarry Smith 
4986356e834SBarry Smith   while (PetscOptionsObject.next) {
4996356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5006356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5016356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5026356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5036356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5046356e834SBarry Smith       } else {
5056356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5066356e834SBarry Smith       }
5076356e834SBarry Smith 
5086356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5096356e834SBarry Smith       case OPTION_HEAD:
5106356e834SBarry Smith         break;
511e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
512e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
513e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
514e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
515e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
516e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
517e26ddf31SBarry Smith         }
518e26ddf31SBarry Smith         break;
5196356e834SBarry Smith       case OPTION_INT:
5207a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5216356e834SBarry Smith         break;
5226356e834SBarry Smith       case OPTION_REAL:
5237a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5246356e834SBarry Smith         break;
5256356e834SBarry Smith       case OPTION_REAL_ARRAY:
5267a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5276356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5287a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5296356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5306356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5316356e834SBarry Smith         }
5326356e834SBarry Smith         break;
5337781c08eSBarry Smith       case OPTION_BOOL:
53471f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5356356e834SBarry Smith         break;
5367781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
537ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5381ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
539ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5401ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5411ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5421ae3d29cSBarry Smith         }
5431ae3d29cSBarry Smith         break;
544a264d7a6SBarry Smith       case OPTION_FLIST:
5451ae3d29cSBarry Smith       case OPTION_ELIST:
5467781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5476356e834SBarry Smith         break;
5481ae3d29cSBarry Smith       case OPTION_STRING:
549475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5508da2146eSBarry Smith         break;
5511ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5521ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5531ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5541ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5551ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5561ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5571ae3d29cSBarry Smith         }
5586356e834SBarry Smith         break;
5596356e834SBarry Smith       }
5606356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5616356e834SBarry Smith     }
562503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
563503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5646356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56571f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
5667781c08eSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
5677781c08eSBarry Smith 
5686356e834SBarry Smith     last                    = PetscOptionsObject.next;
5696356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5706356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5716356e834SBarry Smith   }
5726356e834SBarry Smith   PetscOptionsObject.next = 0;
57353acd3b1SBarry Smith   PetscFunctionReturn(0);
57453acd3b1SBarry Smith }
57553acd3b1SBarry Smith 
57653acd3b1SBarry Smith #undef __FUNCT__
57753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
57853acd3b1SBarry Smith /*@C
57953acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58053acd3b1SBarry Smith 
5813f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58253acd3b1SBarry Smith 
58353acd3b1SBarry Smith    Input Parameters:
58453acd3b1SBarry Smith +  opt - option name
58553acd3b1SBarry Smith .  text - short string that describes the option
58653acd3b1SBarry Smith .  man - manual page with additional information on option
58753acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
58853acd3b1SBarry Smith -  defaultv - the default (current) value
58953acd3b1SBarry Smith 
59053acd3b1SBarry Smith    Output Parameter:
59153acd3b1SBarry Smith +  value - the  value to return
592b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
59353acd3b1SBarry Smith 
59453acd3b1SBarry Smith    Level: beginner
59553acd3b1SBarry Smith 
59653acd3b1SBarry Smith    Concepts: options database
59753acd3b1SBarry Smith 
59853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
59953acd3b1SBarry Smith 
60053acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60153acd3b1SBarry Smith 
60253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
603acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
604acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
607acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
608a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
60953acd3b1SBarry Smith @*/
6107087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61153acd3b1SBarry Smith {
61253acd3b1SBarry Smith   PetscErrorCode ierr;
61353acd3b1SBarry Smith   PetscInt       ntext = 0;
614aa5bb8c0SSatish Balay   PetscInt       tval;
615ace3abfcSBarry Smith   PetscBool      tflg;
61653acd3b1SBarry Smith 
61753acd3b1SBarry Smith   PetscFunctionBegin;
61853acd3b1SBarry Smith   while (list[ntext++]) {
619e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62053acd3b1SBarry Smith   }
621e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62253acd3b1SBarry Smith   ntext -= 3;
623aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
624aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
625aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
626aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
62753acd3b1SBarry Smith   PetscFunctionReturn(0);
62853acd3b1SBarry Smith }
62953acd3b1SBarry Smith 
63053acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63153acd3b1SBarry Smith #undef __FUNCT__
63253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63353acd3b1SBarry Smith /*@C
63453acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63553acd3b1SBarry Smith 
6363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63753acd3b1SBarry Smith 
63853acd3b1SBarry Smith    Input Parameters:
63953acd3b1SBarry Smith +  opt - option name
64053acd3b1SBarry Smith .  text - short string that describes the option
64153acd3b1SBarry Smith .  man - manual page with additional information on option
64253acd3b1SBarry Smith -  defaultv - the default (current) value
64353acd3b1SBarry Smith 
64453acd3b1SBarry Smith    Output Parameter:
64553acd3b1SBarry Smith +  value - the integer value to return
64653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
64753acd3b1SBarry Smith 
64853acd3b1SBarry Smith    Level: beginner
64953acd3b1SBarry Smith 
65053acd3b1SBarry Smith    Concepts: options database^has int
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65353acd3b1SBarry Smith 
65453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
655acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
656acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
659acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
660a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66153acd3b1SBarry Smith @*/
6627087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66353acd3b1SBarry Smith {
66453acd3b1SBarry Smith   PetscErrorCode ierr;
6656356e834SBarry Smith   PetscOptions   amsopt;
66653acd3b1SBarry Smith 
66753acd3b1SBarry Smith   PetscFunctionBegin;
668b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6696356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6706356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
671a297a907SKarl Rupp 
6726356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
673af6d86caSBarry Smith   }
67453acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6762aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
67753acd3b1SBarry Smith   }
67853acd3b1SBarry Smith   PetscFunctionReturn(0);
67953acd3b1SBarry Smith }
68053acd3b1SBarry Smith 
68153acd3b1SBarry Smith #undef __FUNCT__
68253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68353acd3b1SBarry Smith /*@C
68453acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68553acd3b1SBarry Smith 
6863f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith    Input Parameters:
68953acd3b1SBarry Smith +  opt - option name
69053acd3b1SBarry Smith .  text - short string that describes the option
69153acd3b1SBarry Smith .  man - manual page with additional information on option
692bcbf2dc5SJed Brown .  defaultv - the default (current) value
693bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69453acd3b1SBarry Smith 
69553acd3b1SBarry Smith    Output Parameter:
69653acd3b1SBarry Smith +  value - the value to return
69753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69853acd3b1SBarry Smith 
69953acd3b1SBarry Smith    Level: beginner
70053acd3b1SBarry Smith 
70153acd3b1SBarry Smith    Concepts: options database^has int
70253acd3b1SBarry Smith 
70353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70453acd3b1SBarry Smith 
7057fccdfe4SBarry 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).
7067fccdfe4SBarry Smith 
70753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
708acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
709acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
712acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
713a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
71453acd3b1SBarry Smith @*/
7157087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71653acd3b1SBarry Smith {
71753acd3b1SBarry Smith   PetscErrorCode ierr;
718aee2cecaSBarry Smith   PetscOptions   amsopt;
71953acd3b1SBarry Smith 
72053acd3b1SBarry Smith   PetscFunctionBegin;
721b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
722aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
723*77b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
724af6d86caSBarry Smith   }
72553acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
72661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7272aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
72853acd3b1SBarry Smith   }
72953acd3b1SBarry Smith   PetscFunctionReturn(0);
73053acd3b1SBarry Smith }
73153acd3b1SBarry Smith 
73253acd3b1SBarry Smith #undef __FUNCT__
73353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73453acd3b1SBarry Smith /*@C
73553acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
73653acd3b1SBarry Smith 
7373f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
73853acd3b1SBarry Smith 
73953acd3b1SBarry Smith    Input Parameters:
74053acd3b1SBarry Smith +  opt - option name
74153acd3b1SBarry Smith .  text - short string that describes the option
74253acd3b1SBarry Smith .  man - manual page with additional information on option
74353acd3b1SBarry Smith -  defaultv - the default (current) value
74453acd3b1SBarry Smith 
74553acd3b1SBarry Smith    Output Parameter:
74653acd3b1SBarry Smith +  value - the value to return
74753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
74853acd3b1SBarry Smith 
74953acd3b1SBarry Smith    Level: beginner
75053acd3b1SBarry Smith 
75153acd3b1SBarry Smith    Concepts: options database^has int
75253acd3b1SBarry Smith 
75353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75453acd3b1SBarry Smith 
75553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
756acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
757acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
75853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
75953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
760acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
761a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
76253acd3b1SBarry Smith @*/
7637087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76453acd3b1SBarry Smith {
76553acd3b1SBarry Smith   PetscErrorCode ierr;
766538aa990SBarry Smith   PetscOptions   amsopt;
76753acd3b1SBarry Smith 
76853acd3b1SBarry Smith   PetscFunctionBegin;
769b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
770538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
771538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
772a297a907SKarl Rupp 
773538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
774538aa990SBarry Smith   }
77553acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
77661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7772aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
77853acd3b1SBarry Smith   }
77953acd3b1SBarry Smith   PetscFunctionReturn(0);
78053acd3b1SBarry Smith }
78153acd3b1SBarry Smith 
78253acd3b1SBarry Smith #undef __FUNCT__
78353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78453acd3b1SBarry Smith /*@C
78553acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
78653acd3b1SBarry Smith 
7873f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78853acd3b1SBarry Smith 
78953acd3b1SBarry Smith    Input Parameters:
79053acd3b1SBarry Smith +  opt - option name
79153acd3b1SBarry Smith .  text - short string that describes the option
79253acd3b1SBarry Smith .  man - manual page with additional information on option
79353acd3b1SBarry Smith -  defaultv - the default (current) value
79453acd3b1SBarry Smith 
79553acd3b1SBarry Smith    Output Parameter:
79653acd3b1SBarry Smith +  value - the value to return
79753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79853acd3b1SBarry Smith 
79953acd3b1SBarry Smith    Level: beginner
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith    Concepts: options database^has int
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
806acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
807acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
810acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
811a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81253acd3b1SBarry Smith @*/
8137087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81453acd3b1SBarry Smith {
81553acd3b1SBarry Smith   PetscErrorCode ierr;
81653acd3b1SBarry Smith 
81753acd3b1SBarry Smith   PetscFunctionBegin;
81853acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
81953acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82053acd3b1SBarry Smith #else
82153acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82253acd3b1SBarry Smith #endif
82353acd3b1SBarry Smith   PetscFunctionReturn(0);
82453acd3b1SBarry Smith }
82553acd3b1SBarry Smith 
82653acd3b1SBarry Smith #undef __FUNCT__
82753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
82853acd3b1SBarry Smith /*@C
82990d69ab7SBarry 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
83090d69ab7SBarry Smith                       its value is set to false.
83153acd3b1SBarry Smith 
8323f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Input Parameters:
83553acd3b1SBarry Smith +  opt - option name
83653acd3b1SBarry Smith .  text - short string that describes the option
83753acd3b1SBarry Smith -  man - manual page with additional information on option
83853acd3b1SBarry Smith 
83953acd3b1SBarry Smith    Output Parameter:
84053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
84153acd3b1SBarry Smith 
84253acd3b1SBarry Smith    Level: beginner
84353acd3b1SBarry Smith 
84453acd3b1SBarry Smith    Concepts: options database^has int
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
84753acd3b1SBarry Smith 
84853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
849acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
850acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
853acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
854a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
85553acd3b1SBarry Smith @*/
8567087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
85753acd3b1SBarry Smith {
85853acd3b1SBarry Smith   PetscErrorCode ierr;
8591ae3d29cSBarry Smith   PetscOptions   amsopt;
86053acd3b1SBarry Smith 
86153acd3b1SBarry Smith   PetscFunctionBegin;
8621ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8637781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
864ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
865a297a907SKarl Rupp 
866ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8671ae3d29cSBarry Smith   }
86853acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
86961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8702aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
87153acd3b1SBarry Smith   }
87253acd3b1SBarry Smith   PetscFunctionReturn(0);
87353acd3b1SBarry Smith }
87453acd3b1SBarry Smith 
87553acd3b1SBarry Smith #undef __FUNCT__
876a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
87753acd3b1SBarry Smith /*@C
878a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
87953acd3b1SBarry Smith 
8803f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88153acd3b1SBarry Smith 
88253acd3b1SBarry Smith    Input Parameters:
88353acd3b1SBarry Smith +  opt - option name
88453acd3b1SBarry Smith .  text - short string that describes the option
88553acd3b1SBarry Smith .  man - manual page with additional information on option
88653acd3b1SBarry Smith .  list - the possible choices
8873cc1e11dSBarry Smith .  defaultv - the default (current) value
8883cc1e11dSBarry Smith -  len - the length of the character array value
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Output Parameter:
89153acd3b1SBarry Smith +  value - the value to return
89253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith    Level: intermediate
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith    To get a listing of all currently specified options,
90188c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
90253acd3b1SBarry Smith 
903eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
904eabe10d7SBarry Smith 
90553acd3b1SBarry Smith    Concepts: options database^list
90653acd3b1SBarry Smith 
90753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
908acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
90953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
911acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
912a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
91353acd3b1SBarry Smith @*/
914a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91553acd3b1SBarry Smith {
91653acd3b1SBarry Smith   PetscErrorCode ierr;
9173cc1e11dSBarry Smith   PetscOptions   amsopt;
91853acd3b1SBarry Smith 
91953acd3b1SBarry Smith   PetscFunctionBegin;
920b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
921a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
922*77b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
9233cc1e11dSBarry Smith     amsopt->flist = list;
9243cc1e11dSBarry Smith   }
92553acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
92661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
927140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
92853acd3b1SBarry Smith   }
92953acd3b1SBarry Smith   PetscFunctionReturn(0);
93053acd3b1SBarry Smith }
93153acd3b1SBarry Smith 
93253acd3b1SBarry Smith #undef __FUNCT__
93353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
93453acd3b1SBarry Smith /*@C
93553acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
93653acd3b1SBarry Smith 
9373f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93853acd3b1SBarry Smith 
93953acd3b1SBarry Smith    Input Parameters:
94053acd3b1SBarry Smith +  opt - option name
94153acd3b1SBarry Smith .  ltext - short string that describes the option
94253acd3b1SBarry Smith .  man - manual page with additional information on option
943a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
94453acd3b1SBarry Smith .  ntext - number of choices
94553acd3b1SBarry Smith -  defaultv - the default (current) value
94653acd3b1SBarry Smith 
94753acd3b1SBarry Smith    Output Parameter:
94853acd3b1SBarry Smith +  value - the index of the value to return
94953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith    Level: intermediate
95253acd3b1SBarry Smith 
95353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95453acd3b1SBarry Smith 
955a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith    Concepts: options database^list
95853acd3b1SBarry Smith 
95953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
960acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
963acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
964a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
96553acd3b1SBarry Smith @*/
9667087cfbeSBarry 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)
96753acd3b1SBarry Smith {
96853acd3b1SBarry Smith   PetscErrorCode ierr;
96953acd3b1SBarry Smith   PetscInt       i;
9701ae3d29cSBarry Smith   PetscOptions   amsopt;
97153acd3b1SBarry Smith 
97253acd3b1SBarry Smith   PetscFunctionBegin;
9731ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
974d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
975*77b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
9761ae3d29cSBarry Smith     amsopt->list  = list;
9771ae3d29cSBarry Smith     amsopt->nlist = ntext;
9781ae3d29cSBarry Smith   }
97953acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
98153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
98253acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
98353acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
98453acd3b1SBarry Smith     }
9852aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
98653acd3b1SBarry Smith   }
98753acd3b1SBarry Smith   PetscFunctionReturn(0);
98853acd3b1SBarry Smith }
98953acd3b1SBarry Smith 
99053acd3b1SBarry Smith #undef __FUNCT__
991acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
99253acd3b1SBarry Smith /*@C
993acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
994d5649816SBarry Smith        which at most a single value can be true.
99553acd3b1SBarry Smith 
9963f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99753acd3b1SBarry Smith 
99853acd3b1SBarry Smith    Input Parameters:
99953acd3b1SBarry Smith +  opt - option name
100053acd3b1SBarry Smith .  text - short string that describes the option
100153acd3b1SBarry Smith -  man - manual page with additional information on option
100253acd3b1SBarry Smith 
100353acd3b1SBarry Smith    Output Parameter:
100453acd3b1SBarry Smith .  flg - whether that option was set or not
100553acd3b1SBarry Smith 
100653acd3b1SBarry Smith    Level: intermediate
100753acd3b1SBarry Smith 
100853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
100953acd3b1SBarry Smith 
1010acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith     Concepts: options database^logical group
101353acd3b1SBarry Smith 
101453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1015acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1018acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1019a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102053acd3b1SBarry Smith @*/
10217087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
102253acd3b1SBarry Smith {
102353acd3b1SBarry Smith   PetscErrorCode ierr;
10241ae3d29cSBarry Smith   PetscOptions   amsopt;
102553acd3b1SBarry Smith 
102653acd3b1SBarry Smith   PetscFunctionBegin;
10271ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10287781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1029ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1030a297a907SKarl Rupp 
1031ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10321ae3d29cSBarry Smith   }
103368b16fdaSBarry Smith   *flg = PETSC_FALSE;
10340298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
103561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103653acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10372aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
103853acd3b1SBarry Smith   }
103953acd3b1SBarry Smith   PetscFunctionReturn(0);
104053acd3b1SBarry Smith }
104153acd3b1SBarry Smith 
104253acd3b1SBarry Smith #undef __FUNCT__
1043acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
104453acd3b1SBarry Smith /*@C
1045acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1046d5649816SBarry Smith        which at most a single value can be true.
104753acd3b1SBarry Smith 
10483f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104953acd3b1SBarry Smith 
105053acd3b1SBarry Smith    Input Parameters:
105153acd3b1SBarry Smith +  opt - option name
105253acd3b1SBarry Smith .  text - short string that describes the option
105353acd3b1SBarry Smith -  man - manual page with additional information on option
105453acd3b1SBarry Smith 
105553acd3b1SBarry Smith    Output Parameter:
105653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
105753acd3b1SBarry Smith 
105853acd3b1SBarry Smith    Level: intermediate
105953acd3b1SBarry Smith 
106053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106153acd3b1SBarry Smith 
1062acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith     Concepts: options database^logical group
106553acd3b1SBarry Smith 
106653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1067acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
106853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
106953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1070acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1071a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107253acd3b1SBarry Smith @*/
10737087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
107453acd3b1SBarry Smith {
107553acd3b1SBarry Smith   PetscErrorCode ierr;
10761ae3d29cSBarry Smith   PetscOptions   amsopt;
107753acd3b1SBarry Smith 
107853acd3b1SBarry Smith   PetscFunctionBegin;
10791ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10807781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1081ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1082a297a907SKarl Rupp 
1083ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10841ae3d29cSBarry Smith   }
108517326d04SJed Brown   *flg = PETSC_FALSE;
10860298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
108761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10882aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
108953acd3b1SBarry Smith   }
109053acd3b1SBarry Smith   PetscFunctionReturn(0);
109153acd3b1SBarry Smith }
109253acd3b1SBarry Smith 
109353acd3b1SBarry Smith #undef __FUNCT__
1094acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
109553acd3b1SBarry Smith /*@C
1096acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1097d5649816SBarry Smith        which at most a single value can be true.
109853acd3b1SBarry Smith 
10993f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith    Input Parameters:
110253acd3b1SBarry Smith +  opt - option name
110353acd3b1SBarry Smith .  text - short string that describes the option
110453acd3b1SBarry Smith -  man - manual page with additional information on option
110553acd3b1SBarry Smith 
110653acd3b1SBarry Smith    Output Parameter:
110753acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith    Level: intermediate
111053acd3b1SBarry Smith 
111153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111253acd3b1SBarry Smith 
1113acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith     Concepts: options database^logical group
111653acd3b1SBarry Smith 
111753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1118acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
111953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1121acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1122a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
112353acd3b1SBarry Smith @*/
11247087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
112553acd3b1SBarry Smith {
112653acd3b1SBarry Smith   PetscErrorCode ierr;
11271ae3d29cSBarry Smith   PetscOptions   amsopt;
112853acd3b1SBarry Smith 
112953acd3b1SBarry Smith   PetscFunctionBegin;
11301ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11317781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1132ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1133a297a907SKarl Rupp 
1134ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11351ae3d29cSBarry Smith   }
113617326d04SJed Brown   *flg = PETSC_FALSE;
11370298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
113861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11392aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114053acd3b1SBarry Smith   }
114153acd3b1SBarry Smith   PetscFunctionReturn(0);
114253acd3b1SBarry Smith }
114353acd3b1SBarry Smith 
114453acd3b1SBarry Smith #undef __FUNCT__
1145acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
114653acd3b1SBarry Smith /*@C
1147acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
114853acd3b1SBarry Smith 
11493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115053acd3b1SBarry Smith 
115153acd3b1SBarry Smith    Input Parameters:
115253acd3b1SBarry Smith +  opt - option name
115353acd3b1SBarry Smith .  text - short string that describes the option
115453acd3b1SBarry Smith -  man - manual page with additional information on option
115553acd3b1SBarry Smith 
115653acd3b1SBarry Smith    Output Parameter:
115753acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
115853acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
115953acd3b1SBarry Smith 
116053acd3b1SBarry Smith    Level: beginner
116153acd3b1SBarry Smith 
116253acd3b1SBarry Smith    Concepts: options database^logical
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116553acd3b1SBarry Smith 
116653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1167acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1168acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
116953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1171acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1172a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
117353acd3b1SBarry Smith @*/
11747087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
117553acd3b1SBarry Smith {
117653acd3b1SBarry Smith   PetscErrorCode ierr;
1177ace3abfcSBarry Smith   PetscBool      iset;
1178aee2cecaSBarry Smith   PetscOptions   amsopt;
117953acd3b1SBarry Smith 
118053acd3b1SBarry Smith   PetscFunctionBegin;
1181b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
11827781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1183ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1184a297a907SKarl Rupp 
1185ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1186af6d86caSBarry Smith   }
1187acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
118853acd3b1SBarry Smith   if (!iset) {
118953acd3b1SBarry Smith     if (flg) *flg = deflt;
119053acd3b1SBarry Smith   }
119153acd3b1SBarry Smith   if (set) *set = iset;
119261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1193ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
11942aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
119553acd3b1SBarry Smith   }
119653acd3b1SBarry Smith   PetscFunctionReturn(0);
119753acd3b1SBarry Smith }
119853acd3b1SBarry Smith 
119953acd3b1SBarry Smith #undef __FUNCT__
120053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
120153acd3b1SBarry Smith /*@C
120253acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
120353acd3b1SBarry Smith    option in the database. The values must be separated with commas with
120453acd3b1SBarry Smith    no intervening spaces.
120553acd3b1SBarry Smith 
12063f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120753acd3b1SBarry Smith 
120853acd3b1SBarry Smith    Input Parameters:
120953acd3b1SBarry Smith +  opt - the option one is seeking
121053acd3b1SBarry Smith .  text - short string describing option
121153acd3b1SBarry Smith .  man - manual page for option
121253acd3b1SBarry Smith -  nmax - maximum number of values
121353acd3b1SBarry Smith 
121453acd3b1SBarry Smith    Output Parameter:
121553acd3b1SBarry Smith +  value - location to copy values
121653acd3b1SBarry Smith .  nmax - actual number of values found
121753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
121853acd3b1SBarry Smith 
121953acd3b1SBarry Smith    Level: beginner
122053acd3b1SBarry Smith 
122153acd3b1SBarry Smith    Notes:
122253acd3b1SBarry Smith    The user should pass in an array of doubles
122353acd3b1SBarry Smith 
122453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122553acd3b1SBarry Smith 
122653acd3b1SBarry Smith    Concepts: options database^array of strings
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1229acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1232acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1233a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
123453acd3b1SBarry Smith @*/
12357087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
123653acd3b1SBarry Smith {
123753acd3b1SBarry Smith   PetscErrorCode ierr;
123853acd3b1SBarry Smith   PetscInt       i;
1239e26ddf31SBarry Smith   PetscOptions   amsopt;
124053acd3b1SBarry Smith 
124153acd3b1SBarry Smith   PetscFunctionBegin;
1242b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1243e26ddf31SBarry Smith     PetscReal *vals;
1244e26ddf31SBarry Smith 
1245e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1246e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1247e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1248e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1249e26ddf31SBarry Smith     amsopt->arraylength = *n;
1250e26ddf31SBarry Smith   }
125153acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
125261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1253a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
125453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1255a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
125653acd3b1SBarry Smith     }
12572aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
125853acd3b1SBarry Smith   }
125953acd3b1SBarry Smith   PetscFunctionReturn(0);
126053acd3b1SBarry Smith }
126153acd3b1SBarry Smith 
126253acd3b1SBarry Smith 
126353acd3b1SBarry Smith #undef __FUNCT__
126453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
126553acd3b1SBarry Smith /*@C
126653acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1267b32a342fSShri Abhyankar    option in the database.
126853acd3b1SBarry Smith 
12693f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127053acd3b1SBarry Smith 
127153acd3b1SBarry Smith    Input Parameters:
127253acd3b1SBarry Smith +  opt - the option one is seeking
127353acd3b1SBarry Smith .  text - short string describing option
127453acd3b1SBarry Smith .  man - manual page for option
1275f8a50e2bSBarry Smith -  n - maximum number of values
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Output Parameter:
127853acd3b1SBarry Smith +  value - location to copy values
1279f8a50e2bSBarry Smith .  n - actual number of values found
128053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith    Level: beginner
128353acd3b1SBarry Smith 
128453acd3b1SBarry Smith    Notes:
1285b32a342fSShri Abhyankar    The array can be passed as
1286b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12870fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12880fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12890fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1290b32a342fSShri Abhyankar 
1291b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
129253acd3b1SBarry Smith 
129353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129453acd3b1SBarry Smith 
1295b32a342fSShri Abhyankar    Concepts: options database^array of ints
129653acd3b1SBarry Smith 
129753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1298acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1301acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1302a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
130353acd3b1SBarry Smith @*/
13047087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
130553acd3b1SBarry Smith {
130653acd3b1SBarry Smith   PetscErrorCode ierr;
130753acd3b1SBarry Smith   PetscInt       i;
1308e26ddf31SBarry Smith   PetscOptions   amsopt;
130953acd3b1SBarry Smith 
131053acd3b1SBarry Smith   PetscFunctionBegin;
1311b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1312e26ddf31SBarry Smith     PetscInt *vals;
1313e26ddf31SBarry Smith 
1314e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1315e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1316e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1317e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1318e26ddf31SBarry Smith     amsopt->arraylength = *n;
1319e26ddf31SBarry Smith   }
132053acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
132161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
132353acd3b1SBarry Smith     for (i=1; i<*n; i++) {
132453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
132553acd3b1SBarry Smith     }
13262aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
132753acd3b1SBarry Smith   }
132853acd3b1SBarry Smith   PetscFunctionReturn(0);
132953acd3b1SBarry Smith }
133053acd3b1SBarry Smith 
133153acd3b1SBarry Smith #undef __FUNCT__
133253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
133353acd3b1SBarry Smith /*@C
133453acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
133553acd3b1SBarry Smith    option in the database. The values must be separated with commas with
133653acd3b1SBarry Smith    no intervening spaces.
133753acd3b1SBarry Smith 
13383f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
133953acd3b1SBarry Smith 
134053acd3b1SBarry Smith    Input Parameters:
134153acd3b1SBarry Smith +  opt - the option one is seeking
134253acd3b1SBarry Smith .  text - short string describing option
134353acd3b1SBarry Smith .  man - manual page for option
134453acd3b1SBarry Smith -  nmax - maximum number of strings
134553acd3b1SBarry Smith 
134653acd3b1SBarry Smith    Output Parameter:
134753acd3b1SBarry Smith +  value - location to copy strings
134853acd3b1SBarry Smith .  nmax - actual number of strings found
134953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135053acd3b1SBarry Smith 
135153acd3b1SBarry Smith    Level: beginner
135253acd3b1SBarry Smith 
135353acd3b1SBarry Smith    Notes:
135453acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
135553acd3b1SBarry Smith    strings returned by this function.
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
135853acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136153acd3b1SBarry Smith 
136253acd3b1SBarry Smith    Concepts: options database^array of strings
136353acd3b1SBarry Smith 
136453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1365acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
136653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
136753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1368acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1369a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
137053acd3b1SBarry Smith @*/
13717087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
137253acd3b1SBarry Smith {
137353acd3b1SBarry Smith   PetscErrorCode ierr;
13741ae3d29cSBarry Smith   PetscOptions   amsopt;
137553acd3b1SBarry Smith 
137653acd3b1SBarry Smith   PetscFunctionBegin;
13771ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13781ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13791ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1380a297a907SKarl Rupp 
13811ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13821ae3d29cSBarry Smith   }
138353acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
138461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13852aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
138653acd3b1SBarry Smith   }
138753acd3b1SBarry Smith   PetscFunctionReturn(0);
138853acd3b1SBarry Smith }
138953acd3b1SBarry Smith 
1390e2446a98SMatthew Knepley #undef __FUNCT__
1391acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1392e2446a98SMatthew Knepley /*@C
1393acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1394e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1395e2446a98SMatthew Knepley    no intervening spaces.
1396e2446a98SMatthew Knepley 
13973f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1398e2446a98SMatthew Knepley 
1399e2446a98SMatthew Knepley    Input Parameters:
1400e2446a98SMatthew Knepley +  opt - the option one is seeking
1401e2446a98SMatthew Knepley .  text - short string describing option
1402e2446a98SMatthew Knepley .  man - manual page for option
1403e2446a98SMatthew Knepley -  nmax - maximum number of values
1404e2446a98SMatthew Knepley 
1405e2446a98SMatthew Knepley    Output Parameter:
1406e2446a98SMatthew Knepley +  value - location to copy values
1407e2446a98SMatthew Knepley .  nmax - actual number of values found
1408e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1409e2446a98SMatthew Knepley 
1410e2446a98SMatthew Knepley    Level: beginner
1411e2446a98SMatthew Knepley 
1412e2446a98SMatthew Knepley    Notes:
1413e2446a98SMatthew Knepley    The user should pass in an array of doubles
1414e2446a98SMatthew Knepley 
1415e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1416e2446a98SMatthew Knepley 
1417e2446a98SMatthew Knepley    Concepts: options database^array of strings
1418e2446a98SMatthew Knepley 
1419e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1420acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1421e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1422e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1423acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1424a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1425e2446a98SMatthew Knepley @*/
14267087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1427e2446a98SMatthew Knepley {
1428e2446a98SMatthew Knepley   PetscErrorCode ierr;
1429e2446a98SMatthew Knepley   PetscInt       i;
14301ae3d29cSBarry Smith   PetscOptions   amsopt;
1431e2446a98SMatthew Knepley 
1432e2446a98SMatthew Knepley   PetscFunctionBegin;
14331ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1434ace3abfcSBarry Smith     PetscBool *vals;
14351ae3d29cSBarry Smith 
14367781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1437ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1438ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14391ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14401ae3d29cSBarry Smith     amsopt->arraylength = *n;
14411ae3d29cSBarry Smith   }
1442acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1443e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1444e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1445e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1446e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1447e2446a98SMatthew Knepley     }
14482aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1449e2446a98SMatthew Knepley   }
1450e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1451e2446a98SMatthew Knepley }
1452e2446a98SMatthew Knepley 
14538cc676e6SMatthew G Knepley #undef __FUNCT__
14548cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14558cc676e6SMatthew G Knepley /*@C
14568cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14578cc676e6SMatthew G Knepley 
14588cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14598cc676e6SMatthew G Knepley 
14608cc676e6SMatthew G Knepley    Input Parameters:
14618cc676e6SMatthew G Knepley +  opt - option name
14628cc676e6SMatthew G Knepley .  text - short string that describes the option
14638cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14648cc676e6SMatthew G Knepley 
14658cc676e6SMatthew G Knepley    Output Parameter:
14668cc676e6SMatthew G Knepley +  viewer - the viewer
14678cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14688cc676e6SMatthew G Knepley 
14698cc676e6SMatthew G Knepley    Level: beginner
14708cc676e6SMatthew G Knepley 
14718cc676e6SMatthew G Knepley    Concepts: options database^has int
14728cc676e6SMatthew G Knepley 
14738cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14748cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14758cc676e6SMatthew 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
14768cc676e6SMatthew G Knepley $                                     about the object to standard out
14778cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14788cc676e6SMatthew G Knepley $       draw
14798cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14808cc676e6SMatthew G Knepley 
1481cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14828cc676e6SMatthew G Knepley 
14838cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14848cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14858cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14868cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14878cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14888cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1489a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
14908cc676e6SMatthew G Knepley @*/
1491cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14928cc676e6SMatthew G Knepley {
14938cc676e6SMatthew G Knepley   PetscErrorCode ierr;
14948cc676e6SMatthew G Knepley   PetscOptions   amsopt;
14958cc676e6SMatthew G Knepley 
14968cc676e6SMatthew G Knepley   PetscFunctionBegin;
14978cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
14988cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
1499*77b82169SBarry Smith     ierr = PetscStrallocpy("",(char**) &amsopt->data);CHKERRQ(ierr);
15008cc676e6SMatthew G Knepley   }
1501cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15028cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15038cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15048cc676e6SMatthew G Knepley   }
15058cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15068cc676e6SMatthew G Knepley }
15078cc676e6SMatthew G Knepley 
150853acd3b1SBarry Smith 
150953acd3b1SBarry Smith #undef __FUNCT__
151053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
151153acd3b1SBarry Smith /*@C
1512b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
151353acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
151453acd3b1SBarry Smith 
15153f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
151653acd3b1SBarry Smith 
151753acd3b1SBarry Smith    Input Parameter:
151853acd3b1SBarry Smith .   head - the heading text
151953acd3b1SBarry Smith 
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith    Level: intermediate
152253acd3b1SBarry Smith 
152353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152453acd3b1SBarry Smith 
1525b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
152653acd3b1SBarry Smith 
152753acd3b1SBarry Smith    Concepts: options database^subheading
152853acd3b1SBarry Smith 
152953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1530acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1533acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1534a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
153553acd3b1SBarry Smith @*/
15367087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
153753acd3b1SBarry Smith {
153853acd3b1SBarry Smith   PetscErrorCode ierr;
153953acd3b1SBarry Smith 
154053acd3b1SBarry Smith   PetscFunctionBegin;
154161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
154253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
154353acd3b1SBarry Smith   }
154453acd3b1SBarry Smith   PetscFunctionReturn(0);
154553acd3b1SBarry Smith }
154653acd3b1SBarry Smith 
154753acd3b1SBarry Smith 
154853acd3b1SBarry Smith 
154953acd3b1SBarry Smith 
155053acd3b1SBarry Smith 
155153acd3b1SBarry Smith 
1552