xref: /petsc/src/sys/objects/aoptions.c (revision eabe10d74c99faa247dcbf4c622ace6adf7f5338)
17d0a6c19SBarry Smith 
253acd3b1SBarry Smith /*
33fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
43fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
553acd3b1SBarry Smith 
653acd3b1SBarry Smith */
753acd3b1SBarry Smith 
8afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
9665c2dedSJed Brown #include <petscviewer.h>
1053acd3b1SBarry Smith 
112aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
122aa6d131SJed Brown 
1353acd3b1SBarry Smith /*
1453acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
153fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1653acd3b1SBarry Smith 
1753acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1853acd3b1SBarry Smith */
19f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2053acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
2753acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
3253acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3353acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
346356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3753acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith 
410298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4253acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4361b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4553acd3b1SBarry Smith     }
4661b37b28SSatish Balay   }
4753acd3b1SBarry Smith   PetscFunctionReturn(0);
4853acd3b1SBarry Smith }
4953acd3b1SBarry Smith 
503194b578SJed Brown #undef __FUNCT__
513194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
553194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
563194b578SJed Brown {
573194b578SJed Brown   PetscErrorCode ierr;
583194b578SJed Brown   char           title[256];
593194b578SJed Brown   PetscBool      flg;
603194b578SJed Brown 
613194b578SJed Brown   PetscFunctionBegin;
623194b578SJed Brown   PetscValidHeader(obj,1);
633194b578SJed Brown   PetscOptionsObject.object         = obj;
643194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
65a297a907SKarl Rupp 
663194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
673194b578SJed Brown   if (flg) {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   } else {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   }
723194b578SJed Brown   ierr = PetscOptionsBegin_Private(obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
733194b578SJed Brown   PetscFunctionReturn(0);
743194b578SJed Brown }
753194b578SJed Brown 
7653acd3b1SBarry Smith /*
7753acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7853acd3b1SBarry Smith */
7953acd3b1SBarry Smith #undef __FUNCT__
8053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
81e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8253acd3b1SBarry Smith {
8353acd3b1SBarry Smith   int          ierr;
8453acd3b1SBarry Smith   PetscOptions next;
853be6e4c3SJed Brown   PetscBool    valid;
8653acd3b1SBarry Smith 
8753acd3b1SBarry Smith   PetscFunctionBegin;
883be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
893be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
903be6e4c3SJed Brown 
916e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9253acd3b1SBarry Smith   (*amsopt)->next = 0;
9353acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
946356e834SBarry Smith   (*amsopt)->type = t;
9553acd3b1SBarry Smith   (*amsopt)->data = 0;
9661b37b28SSatish Balay 
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9853acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
996356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10053acd3b1SBarry Smith 
101a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
102a297a907SKarl Rupp   else {
10353acd3b1SBarry Smith     next = PetscOptionsObject.next;
10453acd3b1SBarry Smith     while (next->next) next = next->next;
10553acd3b1SBarry Smith     next->next = *amsopt;
10653acd3b1SBarry Smith   }
10753acd3b1SBarry Smith   PetscFunctionReturn(0);
10853acd3b1SBarry Smith }
10953acd3b1SBarry Smith 
11053acd3b1SBarry Smith #undef __FUNCT__
111aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
112aee2cecaSBarry Smith /*
1133fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith     Collective on MPI_Comm
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith    Input Parameters:
1183fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1193fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1203fc1eb6aSBarry Smith -     str - location to store input
121aee2cecaSBarry Smith 
122aee2cecaSBarry Smith     Bugs:
123aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
124aee2cecaSBarry Smith 
125aee2cecaSBarry Smith */
1263fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
127aee2cecaSBarry Smith {
128330cf3c9SBarry Smith   size_t         i;
129aee2cecaSBarry Smith   char           c;
1303fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
131aee2cecaSBarry Smith   PetscErrorCode ierr;
132aee2cecaSBarry Smith 
133aee2cecaSBarry Smith   PetscFunctionBegin;
134aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
135aee2cecaSBarry Smith   if (!rank) {
136aee2cecaSBarry Smith     c = (char) getchar();
137aee2cecaSBarry Smith     i = 0;
138aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
139aee2cecaSBarry Smith       str[i++] = c;
140aee2cecaSBarry Smith       c = (char)getchar();
141aee2cecaSBarry Smith     }
142aee2cecaSBarry Smith     str[i] = 0;
143aee2cecaSBarry Smith   }
1444dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1453fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
146aee2cecaSBarry Smith   PetscFunctionReturn(0);
147aee2cecaSBarry Smith }
148aee2cecaSBarry Smith 
149aee2cecaSBarry Smith #undef __FUNCT__
150aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
151aee2cecaSBarry Smith /*
1523cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
153aee2cecaSBarry Smith 
154aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
155aee2cecaSBarry Smith 
1567781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1577781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1587781c08eSBarry Smith 
159aee2cecaSBarry Smith     Bugs:
1607781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1613cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
162aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
163aee2cecaSBarry Smith 
1643cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1653cc1e11dSBarry Smith      address space and communicating with the PETSc program
1663cc1e11dSBarry Smith 
167aee2cecaSBarry Smith */
168aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1696356e834SBarry Smith {
1706356e834SBarry Smith   PetscErrorCode ierr;
1716356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1726356e834SBarry Smith   char           str[512];
173a4404d99SBarry Smith   PetscInt       id;
1747781c08eSBarry Smith   PetscBool      bid;
175a4404d99SBarry Smith   PetscReal      ir,*valr;
176330cf3c9SBarry Smith   PetscInt       *vald;
177330cf3c9SBarry Smith   size_t         i;
1786356e834SBarry Smith 
179e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1806356e834SBarry Smith   while (next) {
1816356e834SBarry Smith     switch (next->type) {
1826356e834SBarry Smith     case OPTION_HEAD:
1836356e834SBarry Smith       break;
184e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
185e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
186e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
187e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
188e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
189e26ddf31SBarry Smith         if (i < next->arraylength-1) {
190e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
191e26ddf31SBarry Smith         }
192e26ddf31SBarry Smith       }
193e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
194e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
195e26ddf31SBarry Smith       if (str[0]) {
196e26ddf31SBarry Smith         PetscToken token;
197e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
198e26ddf31SBarry Smith         size_t     len;
199e26ddf31SBarry Smith         char       *value;
200ace3abfcSBarry Smith         PetscBool  foundrange;
201e26ddf31SBarry Smith 
202e26ddf31SBarry Smith         next->set = PETSC_TRUE;
203e26ddf31SBarry Smith         value     = str;
204e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
205e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
206e26ddf31SBarry Smith         while (n < nmax) {
207e26ddf31SBarry Smith           if (!value) break;
208e26ddf31SBarry Smith 
209e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
210e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
211e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
212e26ddf31SBarry Smith           if (value[0] == '-') i=2;
213e26ddf31SBarry Smith           else i=1;
214330cf3c9SBarry Smith           for (;i<len; i++) {
215e26ddf31SBarry Smith             if (value[i] == '-') {
216e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
217e26ddf31SBarry Smith               value[i] = 0;
218cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
219cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
220e32f2f54SBarry Smith               if (end <= start) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, %s-%s cannot have decreasing list",n,value,value+i+1);
221e32f2f54SBarry Smith               if (n + end - start - 1 >= nmax) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, not enough space in left in array (%D) to contain entire range from %D to %D",n,nmax-n,start,end);
222e26ddf31SBarry Smith               for (; start<end; start++) {
223e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
224e26ddf31SBarry Smith               }
225e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
226e26ddf31SBarry Smith               break;
227e26ddf31SBarry Smith             }
228e26ddf31SBarry Smith           }
229e26ddf31SBarry Smith           if (!foundrange) {
230cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
231e26ddf31SBarry Smith             dvalue++;
232e26ddf31SBarry Smith             n++;
233e26ddf31SBarry Smith           }
234e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
235e26ddf31SBarry Smith         }
2368c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
237e26ddf31SBarry Smith       }
238e26ddf31SBarry Smith       break;
239e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
240e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
241e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
242e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
243e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
244e26ddf31SBarry Smith         if (i < next->arraylength-1) {
245e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
246e26ddf31SBarry Smith         }
247e26ddf31SBarry Smith       }
248e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
249e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
250e26ddf31SBarry Smith       if (str[0]) {
251e26ddf31SBarry Smith         PetscToken token;
252e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
253e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
254e26ddf31SBarry Smith         char       *value;
255e26ddf31SBarry Smith 
256e26ddf31SBarry Smith         next->set = PETSC_TRUE;
257e26ddf31SBarry Smith         value     = str;
258e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
259e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
260e26ddf31SBarry Smith         while (n < nmax) {
261e26ddf31SBarry Smith           if (!value) break;
262cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
263e26ddf31SBarry Smith           dvalue++;
264e26ddf31SBarry Smith           n++;
265e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
266e26ddf31SBarry Smith         }
2678c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
268e26ddf31SBarry Smith       }
269e26ddf31SBarry Smith       break;
2706356e834SBarry Smith     case OPTION_INT:
271e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%d>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(int*)next->data,next->text,next->man);CHKERRQ(ierr);
2723fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2733fc1eb6aSBarry Smith       if (str[0]) {
274c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
275c272547aSJed Brown         sscanf(str,"%lld",&id);
276c272547aSJed Brown #else
277aee2cecaSBarry Smith         sscanf(str,"%d",&id);
278c272547aSJed Brown #endif
279aee2cecaSBarry Smith         next->set = PETSC_TRUE;
280a297a907SKarl Rupp 
281aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
282aee2cecaSBarry Smith       }
283aee2cecaSBarry Smith       break;
284aee2cecaSBarry Smith     case OPTION_REAL:
285e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%g>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(double*)next->data,next->text,next->man);CHKERRQ(ierr);
2863fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2873fc1eb6aSBarry Smith       if (str[0]) {
288ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
289a4404d99SBarry Smith         sscanf(str,"%e",&ir);
290ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
291aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
292ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
293d9822059SBarry Smith         ir = strtoflt128(str,0);
294d9822059SBarry Smith #else
295513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
296a4404d99SBarry Smith #endif
297aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
298aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
299aee2cecaSBarry Smith       }
300aee2cecaSBarry Smith       break;
3017781c08eSBarry Smith     case OPTION_BOOL:
3027781c08eSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3037781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3047781c08eSBarry Smith       if (str[0]) {
3057781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3067781c08eSBarry Smith         next->set = PETSC_TRUE;
3077781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3087781c08eSBarry Smith       }
3097781c08eSBarry Smith       break;
310aee2cecaSBarry Smith     case OPTION_STRING:
311e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,(char*)next->data,next->text,next->man);CHKERRQ(ierr);
3123fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3133fc1eb6aSBarry Smith       if (str[0]) {
314aee2cecaSBarry Smith         next->set = PETSC_TRUE;
3157781c08eSBarry Smith         next->data = (void*)strdup(str);
3166356e834SBarry Smith       }
3176356e834SBarry Smith       break;
3183cc1e11dSBarry Smith     case OPTION_LIST:
319140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3203cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3213cc1e11dSBarry Smith       if (str[0]) {
3223cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3233cc1e11dSBarry Smith         next->set = PETSC_TRUE;
3247781c08eSBarry Smith         next->data = (void*)strdup(str);
3253cc1e11dSBarry Smith       }
3263cc1e11dSBarry Smith       break;
327b432afdaSMatthew Knepley     default:
328b432afdaSMatthew Knepley       break;
3296356e834SBarry Smith     }
3306356e834SBarry Smith     next = next->next;
3316356e834SBarry Smith   }
3326356e834SBarry Smith   PetscFunctionReturn(0);
3336356e834SBarry Smith }
3346356e834SBarry Smith 
335e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
336e04113cfSBarry Smith #include <petscviewersaws.h>
337d5649816SBarry Smith 
338d5649816SBarry Smith static int count = 0;
339d5649816SBarry Smith 
340b3506946SBarry Smith #undef __FUNCT__
341e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
342e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
343d5649816SBarry Smith {
3442657e9d9SBarry Smith   PetscFunctionBegin;
345d5649816SBarry Smith   PetscFunctionReturn(0);
346d5649816SBarry Smith }
347d5649816SBarry Smith 
348d5649816SBarry Smith #undef __FUNCT__
3497781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
350b3506946SBarry Smith /*
3517781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
352b3506946SBarry Smith 
353b3506946SBarry Smith     Bugs:
354b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
355b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
356b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
357b3506946SBarry Smith 
358b3506946SBarry Smith 
359b3506946SBarry Smith */
360475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
361b3506946SBarry Smith {
362b3506946SBarry Smith   PetscErrorCode ierr;
363b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
364d5649816SBarry Smith   static int     mancount = 0;
365b3506946SBarry Smith   char           options[16];
3667aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
3677781c08eSBarry Smith   char           manname[16],setname[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 
3757781c08eSBarry 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) {
3827781c08eSBarry Smith     sprintf(setname,"set_%d",mancount);
3837781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",setname);CHKERRQ(ierr);
3842657e9d9SBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->set,1,SAWs_WRITE,SAWs_INT));
3857781c08eSBarry Smith     sprintf(manname,"man_%d",mancount);
3862657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
3877781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
3887781c08eSBarry Smith     sprintf(textname,"text_%d",mancount++);
3897781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
3907781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
3919f32e415SBarry Smith 
392b3506946SBarry Smith     switch (next->type) {
393b3506946SBarry Smith     case OPTION_HEAD:
394b3506946SBarry Smith       break;
395b3506946SBarry Smith     case OPTION_INT_ARRAY:
3967781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
3972657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
398b3506946SBarry Smith       break;
399b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4007781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4012657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
402b3506946SBarry Smith       break;
403b3506946SBarry Smith     case OPTION_INT:
4047781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4052657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
406b3506946SBarry Smith       break;
407b3506946SBarry Smith     case OPTION_REAL:
4087781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4092657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
410b3506946SBarry Smith       break;
4117781c08eSBarry Smith     case OPTION_BOOL:
4127781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4132657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4141ae3d29cSBarry Smith       break;
4157781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4167781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4172657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
41871f08665SBarry Smith       break;
419b3506946SBarry Smith     case OPTION_STRING:
4207781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4217781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4221ae3d29cSBarry Smith       break;
4231ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4247781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4252657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
426b3506946SBarry Smith       break;
427b3506946SBarry Smith     case OPTION_LIST:
4287781c08eSBarry Smith       {/*PetscInt ntext;*/
4297781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4307781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4317781c08eSBarry Smith       /*ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
4327781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4337781c08eSBarry Smith        PetscStackCallSAWs(SAWs_Register,(dir,next->edata,ntext-1,SAWs_READ,SAWs_STRING));*/
4341ae3d29cSBarry Smith       break;}
4351ae3d29cSBarry Smith     case OPTION_ELIST:
4367781c08eSBarry Smith       {/*PetscInt ntext = next->nlist; */
4377781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4387781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4397781c08eSBarry Smith       /*ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4401ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4412657e9d9SBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->text);CHKERRQ(ierr);
4427781c08eSBarry Smith        PetscStackCallSAWs(SAWs_Register,(dir,next->edata,ntext,SAWs_READ,SAWs_STRING));*/
44371f08665SBarry Smith       break;}
444b3506946SBarry Smith     default:
445b3506946SBarry Smith       break;
446b3506946SBarry Smith     }
447b3506946SBarry Smith     next = next->next;
448b3506946SBarry Smith   }
449b3506946SBarry Smith 
450b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4517aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
452b3506946SBarry Smith 
453b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
454b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
455b3506946SBarry Smith 
4569a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
457b3506946SBarry Smith   PetscFunctionReturn(0);
458b3506946SBarry Smith }
459b3506946SBarry Smith #endif
460b3506946SBarry Smith 
4616356e834SBarry Smith #undef __FUNCT__
46253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
46353acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
46453acd3b1SBarry Smith {
46553acd3b1SBarry Smith   PetscErrorCode ierr;
4666356e834SBarry Smith   PetscOptions   last;
4676356e834SBarry Smith   char           option[256],value[1024],tmp[32];
468330cf3c9SBarry Smith   size_t         j;
46953acd3b1SBarry Smith 
47053acd3b1SBarry Smith   PetscFunctionBegin;
471aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
472b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
4737781c08eSBarry Smith #if defined(PETSC_HAVE_SAWS)
474475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
475b3506946SBarry Smith #else
47671f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
477b3506946SBarry Smith #endif
478aee2cecaSBarry Smith     }
479aee2cecaSBarry Smith   }
4806356e834SBarry Smith 
481c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
482c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4836356e834SBarry Smith 
4846356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4856356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
48661b37b28SSatish Balay   /* reset alreadyprinted flag */
48761b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4883194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4890298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4906356e834SBarry Smith 
4916356e834SBarry Smith   while (PetscOptionsObject.next) {
4926356e834SBarry Smith     if (PetscOptionsObject.next->set) {
4936356e834SBarry Smith       if (PetscOptionsObject.prefix) {
4946356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
4956356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
4966356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
4976356e834SBarry Smith       } else {
4986356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
4996356e834SBarry Smith       }
5006356e834SBarry Smith 
5016356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5026356e834SBarry Smith       case OPTION_HEAD:
5036356e834SBarry Smith         break;
504e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
505e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
506e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
507e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
508e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
509e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
510e26ddf31SBarry Smith         }
511e26ddf31SBarry Smith         break;
5126356e834SBarry Smith       case OPTION_INT:
5137a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5146356e834SBarry Smith         break;
5156356e834SBarry Smith       case OPTION_REAL:
5167a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5176356e834SBarry Smith         break;
5186356e834SBarry Smith       case OPTION_REAL_ARRAY:
5197a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5206356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5217a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5226356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5236356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5246356e834SBarry Smith         }
5256356e834SBarry Smith         break;
5267781c08eSBarry Smith       case OPTION_BOOL:
52771f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5286356e834SBarry Smith         break;
5297781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
530ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5311ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
532ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5331ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5341ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5351ae3d29cSBarry Smith         }
5361ae3d29cSBarry Smith         break;
5376356e834SBarry Smith       case OPTION_LIST:
5381ae3d29cSBarry Smith       case OPTION_ELIST:
5397781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5406356e834SBarry Smith         break;
5411ae3d29cSBarry Smith       case OPTION_STRING:
542475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5431ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5441ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5451ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5461ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5471ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5481ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5491ae3d29cSBarry Smith         }
5506356e834SBarry Smith         break;
5516356e834SBarry Smith       }
5526356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5536356e834SBarry Smith     }
554503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
555503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5566356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
55771f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
558a297a907SKarl Rupp 
5597781c08eSBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_LIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
5607781c08eSBarry Smith       free(PetscOptionsObject.next->data);
5617781c08eSBarry Smith     } else {
5627781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
5637781c08eSBarry Smith     }
5647781c08eSBarry Smith 
5656356e834SBarry Smith     last                    = PetscOptionsObject.next;
5666356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5676356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5686356e834SBarry Smith   }
5696356e834SBarry Smith   PetscOptionsObject.next = 0;
57053acd3b1SBarry Smith   PetscFunctionReturn(0);
57153acd3b1SBarry Smith }
57253acd3b1SBarry Smith 
57353acd3b1SBarry Smith #undef __FUNCT__
57453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
57553acd3b1SBarry Smith /*@C
57653acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
57753acd3b1SBarry Smith 
5783f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
57953acd3b1SBarry Smith 
58053acd3b1SBarry Smith    Input Parameters:
58153acd3b1SBarry Smith +  opt - option name
58253acd3b1SBarry Smith .  text - short string that describes the option
58353acd3b1SBarry Smith .  man - manual page with additional information on option
58453acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
58553acd3b1SBarry Smith -  defaultv - the default (current) value
58653acd3b1SBarry Smith 
58753acd3b1SBarry Smith    Output Parameter:
58853acd3b1SBarry Smith +  value - the  value to return
589b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
59053acd3b1SBarry Smith 
59153acd3b1SBarry Smith    Level: beginner
59253acd3b1SBarry Smith 
59353acd3b1SBarry Smith    Concepts: options database
59453acd3b1SBarry Smith 
59553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
59653acd3b1SBarry Smith 
59753acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
59853acd3b1SBarry Smith 
59953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
600acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
601acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
604acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
60553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
60653acd3b1SBarry Smith @*/
6077087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
60853acd3b1SBarry Smith {
60953acd3b1SBarry Smith   PetscErrorCode ierr;
61053acd3b1SBarry Smith   PetscInt       ntext = 0;
611aa5bb8c0SSatish Balay   PetscInt       tval;
612ace3abfcSBarry Smith   PetscBool      tflg;
61353acd3b1SBarry Smith 
61453acd3b1SBarry Smith   PetscFunctionBegin;
61553acd3b1SBarry Smith   while (list[ntext++]) {
616e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
61753acd3b1SBarry Smith   }
618e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
61953acd3b1SBarry Smith   ntext -= 3;
620aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
621aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
622aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
623aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
62453acd3b1SBarry Smith   PetscFunctionReturn(0);
62553acd3b1SBarry Smith }
62653acd3b1SBarry Smith 
62753acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
62853acd3b1SBarry Smith #undef __FUNCT__
62953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63053acd3b1SBarry Smith /*@C
63153acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63253acd3b1SBarry Smith 
6333f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63453acd3b1SBarry Smith 
63553acd3b1SBarry Smith    Input Parameters:
63653acd3b1SBarry Smith +  opt - option name
63753acd3b1SBarry Smith .  text - short string that describes the option
63853acd3b1SBarry Smith .  man - manual page with additional information on option
63953acd3b1SBarry Smith -  defaultv - the default (current) value
64053acd3b1SBarry Smith 
64153acd3b1SBarry Smith    Output Parameter:
64253acd3b1SBarry Smith +  value - the integer value to return
64353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
64453acd3b1SBarry Smith 
64553acd3b1SBarry Smith    Level: beginner
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith    Concepts: options database^has int
64853acd3b1SBarry Smith 
64953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
652acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
653acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
656acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
65753acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
65853acd3b1SBarry Smith @*/
6597087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66053acd3b1SBarry Smith {
66153acd3b1SBarry Smith   PetscErrorCode ierr;
6626356e834SBarry Smith   PetscOptions   amsopt;
66353acd3b1SBarry Smith 
66453acd3b1SBarry Smith   PetscFunctionBegin;
665b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6666356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6676356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
668a297a907SKarl Rupp 
6696356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
670af6d86caSBarry Smith   }
67153acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6732aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
67453acd3b1SBarry Smith   }
67553acd3b1SBarry Smith   PetscFunctionReturn(0);
67653acd3b1SBarry Smith }
67753acd3b1SBarry Smith 
67853acd3b1SBarry Smith #undef __FUNCT__
67953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68053acd3b1SBarry Smith /*@C
68153acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68253acd3b1SBarry Smith 
6833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68453acd3b1SBarry Smith 
68553acd3b1SBarry Smith    Input Parameters:
68653acd3b1SBarry Smith +  opt - option name
68753acd3b1SBarry Smith .  text - short string that describes the option
68853acd3b1SBarry Smith .  man - manual page with additional information on option
689bcbf2dc5SJed Brown .  defaultv - the default (current) value
690bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith    Output Parameter:
69353acd3b1SBarry Smith +  value - the value to return
69453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith    Level: beginner
69753acd3b1SBarry Smith 
69853acd3b1SBarry Smith    Concepts: options database^has int
69953acd3b1SBarry Smith 
70053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70153acd3b1SBarry Smith 
7027fccdfe4SBarry 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).
7037fccdfe4SBarry Smith 
70453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
705acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
706acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
70753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
70853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
709acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71153acd3b1SBarry Smith @*/
7127087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71353acd3b1SBarry Smith {
71453acd3b1SBarry Smith   PetscErrorCode ierr;
715aee2cecaSBarry Smith   PetscOptions   amsopt;
71653acd3b1SBarry Smith 
71753acd3b1SBarry Smith   PetscFunctionBegin;
718b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
719aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
7207781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
721af6d86caSBarry Smith   }
72253acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
72361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7242aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
72553acd3b1SBarry Smith   }
72653acd3b1SBarry Smith   PetscFunctionReturn(0);
72753acd3b1SBarry Smith }
72853acd3b1SBarry Smith 
72953acd3b1SBarry Smith #undef __FUNCT__
73053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73153acd3b1SBarry Smith /*@C
73253acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
73353acd3b1SBarry Smith 
7343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
73553acd3b1SBarry Smith 
73653acd3b1SBarry Smith    Input Parameters:
73753acd3b1SBarry Smith +  opt - option name
73853acd3b1SBarry Smith .  text - short string that describes the option
73953acd3b1SBarry Smith .  man - manual page with additional information on option
74053acd3b1SBarry Smith -  defaultv - the default (current) value
74153acd3b1SBarry Smith 
74253acd3b1SBarry Smith    Output Parameter:
74353acd3b1SBarry Smith +  value - the value to return
74453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
74553acd3b1SBarry Smith 
74653acd3b1SBarry Smith    Level: beginner
74753acd3b1SBarry Smith 
74853acd3b1SBarry Smith    Concepts: options database^has int
74953acd3b1SBarry Smith 
75053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75153acd3b1SBarry Smith 
75253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
753acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
754acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
75553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
75653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
757acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
75853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
75953acd3b1SBarry Smith @*/
7607087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76153acd3b1SBarry Smith {
76253acd3b1SBarry Smith   PetscErrorCode ierr;
763538aa990SBarry Smith   PetscOptions   amsopt;
76453acd3b1SBarry Smith 
76553acd3b1SBarry Smith   PetscFunctionBegin;
766b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
767538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
768538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
769a297a907SKarl Rupp 
770538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
771538aa990SBarry Smith   }
77253acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
77361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7742aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
77553acd3b1SBarry Smith   }
77653acd3b1SBarry Smith   PetscFunctionReturn(0);
77753acd3b1SBarry Smith }
77853acd3b1SBarry Smith 
77953acd3b1SBarry Smith #undef __FUNCT__
78053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78153acd3b1SBarry Smith /*@C
78253acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
78353acd3b1SBarry Smith 
7843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
78553acd3b1SBarry Smith 
78653acd3b1SBarry Smith    Input Parameters:
78753acd3b1SBarry Smith +  opt - option name
78853acd3b1SBarry Smith .  text - short string that describes the option
78953acd3b1SBarry Smith .  man - manual page with additional information on option
79053acd3b1SBarry Smith -  defaultv - the default (current) value
79153acd3b1SBarry Smith 
79253acd3b1SBarry Smith    Output Parameter:
79353acd3b1SBarry Smith +  value - the value to return
79453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith    Level: beginner
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith    Concepts: options database^has int
79953acd3b1SBarry Smith 
80053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80153acd3b1SBarry Smith 
80253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
803acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
804acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
807acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
80853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
80953acd3b1SBarry Smith @*/
8107087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81153acd3b1SBarry Smith {
81253acd3b1SBarry Smith   PetscErrorCode ierr;
81353acd3b1SBarry Smith 
81453acd3b1SBarry Smith   PetscFunctionBegin;
81553acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
81653acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
81753acd3b1SBarry Smith #else
81853acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
81953acd3b1SBarry Smith #endif
82053acd3b1SBarry Smith   PetscFunctionReturn(0);
82153acd3b1SBarry Smith }
82253acd3b1SBarry Smith 
82353acd3b1SBarry Smith #undef __FUNCT__
82453acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
82553acd3b1SBarry Smith /*@C
82690d69ab7SBarry 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
82790d69ab7SBarry Smith                       its value is set to false.
82853acd3b1SBarry Smith 
8293f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83053acd3b1SBarry Smith 
83153acd3b1SBarry Smith    Input Parameters:
83253acd3b1SBarry Smith +  opt - option name
83353acd3b1SBarry Smith .  text - short string that describes the option
83453acd3b1SBarry Smith -  man - manual page with additional information on option
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Output Parameter:
83753acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
83853acd3b1SBarry Smith 
83953acd3b1SBarry Smith    Level: beginner
84053acd3b1SBarry Smith 
84153acd3b1SBarry Smith    Concepts: options database^has int
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
84453acd3b1SBarry Smith 
84553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
846acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
847acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
84853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
850acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
85153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
85253acd3b1SBarry Smith @*/
8537087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
85453acd3b1SBarry Smith {
85553acd3b1SBarry Smith   PetscErrorCode ierr;
8561ae3d29cSBarry Smith   PetscOptions   amsopt;
85753acd3b1SBarry Smith 
85853acd3b1SBarry Smith   PetscFunctionBegin;
8591ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8607781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
861ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
862a297a907SKarl Rupp 
863ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8641ae3d29cSBarry Smith   }
86553acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
86661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8672aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
86853acd3b1SBarry Smith   }
86953acd3b1SBarry Smith   PetscFunctionReturn(0);
87053acd3b1SBarry Smith }
87153acd3b1SBarry Smith 
87253acd3b1SBarry Smith #undef __FUNCT__
87353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
87453acd3b1SBarry Smith /*@C
87553acd3b1SBarry Smith      PetscOptionsList - Puts a list of option values that a single one may be selected from
87653acd3b1SBarry Smith 
8773f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
87853acd3b1SBarry Smith 
87953acd3b1SBarry Smith    Input Parameters:
88053acd3b1SBarry Smith +  opt - option name
88153acd3b1SBarry Smith .  text - short string that describes the option
88253acd3b1SBarry Smith .  man - manual page with additional information on option
88353acd3b1SBarry Smith .  list - the possible choices
8843cc1e11dSBarry Smith .  defaultv - the default (current) value
8853cc1e11dSBarry Smith -  len - the length of the character array value
88653acd3b1SBarry Smith 
88753acd3b1SBarry Smith    Output Parameter:
88853acd3b1SBarry Smith +  value - the value to return
88953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89053acd3b1SBarry Smith 
89153acd3b1SBarry Smith    Level: intermediate
89253acd3b1SBarry Smith 
89353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
89653acd3b1SBarry Smith 
89753acd3b1SBarry Smith    To get a listing of all currently specified options,
89888c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
89953acd3b1SBarry Smith 
900*eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
901*eabe10d7SBarry Smith 
90253acd3b1SBarry Smith    Concepts: options database^list
90353acd3b1SBarry Smith 
90453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
905acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
90653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
908acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
909171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
91053acd3b1SBarry Smith @*/
911140e18c1SBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91253acd3b1SBarry Smith {
91353acd3b1SBarry Smith   PetscErrorCode ierr;
9143cc1e11dSBarry Smith   PetscOptions   amsopt;
91553acd3b1SBarry Smith 
91653acd3b1SBarry Smith   PetscFunctionBegin;
917b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9183cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
9197781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9203cc1e11dSBarry Smith     amsopt->flist = list;
9213cc1e11dSBarry Smith   }
92253acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
92361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
924140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
92553acd3b1SBarry Smith   }
92653acd3b1SBarry Smith   PetscFunctionReturn(0);
92753acd3b1SBarry Smith }
92853acd3b1SBarry Smith 
92953acd3b1SBarry Smith #undef __FUNCT__
93053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
93153acd3b1SBarry Smith /*@C
93253acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
93353acd3b1SBarry Smith 
9343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
93553acd3b1SBarry Smith 
93653acd3b1SBarry Smith    Input Parameters:
93753acd3b1SBarry Smith +  opt - option name
93853acd3b1SBarry Smith .  ltext - short string that describes the option
93953acd3b1SBarry Smith .  man - manual page with additional information on option
94053acd3b1SBarry Smith .  list - the possible choices
94153acd3b1SBarry Smith .  ntext - number of choices
94253acd3b1SBarry Smith -  defaultv - the default (current) value
94353acd3b1SBarry Smith 
94453acd3b1SBarry Smith    Output Parameter:
94553acd3b1SBarry Smith +  value - the index of the value to return
94653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94753acd3b1SBarry Smith 
94853acd3b1SBarry Smith    Level: intermediate
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95153acd3b1SBarry Smith 
952140e18c1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFunctionList()
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith    Concepts: options database^list
95553acd3b1SBarry Smith 
95653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
957acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
95853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
95953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
960acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
961171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
96253acd3b1SBarry Smith @*/
9637087cfbeSBarry 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)
96453acd3b1SBarry Smith {
96553acd3b1SBarry Smith   PetscErrorCode ierr;
96653acd3b1SBarry Smith   PetscInt       i;
9671ae3d29cSBarry Smith   PetscOptions   amsopt;
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith   PetscFunctionBegin;
9701ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
971d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9727781c08eSBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9731ae3d29cSBarry Smith     amsopt->list  = list;
9741ae3d29cSBarry Smith     amsopt->nlist = ntext;
9751ae3d29cSBarry Smith   }
97653acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
97761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
97853acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
97953acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
98053acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
98153acd3b1SBarry Smith     }
9822aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
98353acd3b1SBarry Smith   }
98453acd3b1SBarry Smith   PetscFunctionReturn(0);
98553acd3b1SBarry Smith }
98653acd3b1SBarry Smith 
98753acd3b1SBarry Smith #undef __FUNCT__
988acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
98953acd3b1SBarry Smith /*@C
990acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
991d5649816SBarry Smith        which at most a single value can be true.
99253acd3b1SBarry Smith 
9933f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
99453acd3b1SBarry Smith 
99553acd3b1SBarry Smith    Input Parameters:
99653acd3b1SBarry Smith +  opt - option name
99753acd3b1SBarry Smith .  text - short string that describes the option
99853acd3b1SBarry Smith -  man - manual page with additional information on option
99953acd3b1SBarry Smith 
100053acd3b1SBarry Smith    Output Parameter:
100153acd3b1SBarry Smith .  flg - whether that option was set or not
100253acd3b1SBarry Smith 
100353acd3b1SBarry Smith    Level: intermediate
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
100653acd3b1SBarry Smith 
1007acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith     Concepts: options database^logical group
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1012acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1015acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
101653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
101753acd3b1SBarry Smith @*/
10187087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
101953acd3b1SBarry Smith {
102053acd3b1SBarry Smith   PetscErrorCode ierr;
10211ae3d29cSBarry Smith   PetscOptions   amsopt;
102253acd3b1SBarry Smith 
102353acd3b1SBarry Smith   PetscFunctionBegin;
10241ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10257781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1026ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1027a297a907SKarl Rupp 
1028ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10291ae3d29cSBarry Smith   }
103068b16fdaSBarry Smith   *flg = PETSC_FALSE;
10310298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
103261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103353acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10342aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
103553acd3b1SBarry Smith   }
103653acd3b1SBarry Smith   PetscFunctionReturn(0);
103753acd3b1SBarry Smith }
103853acd3b1SBarry Smith 
103953acd3b1SBarry Smith #undef __FUNCT__
1040acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
104153acd3b1SBarry Smith /*@C
1042acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1043d5649816SBarry Smith        which at most a single value can be true.
104453acd3b1SBarry Smith 
10453f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
104653acd3b1SBarry Smith 
104753acd3b1SBarry Smith    Input Parameters:
104853acd3b1SBarry Smith +  opt - option name
104953acd3b1SBarry Smith .  text - short string that describes the option
105053acd3b1SBarry Smith -  man - manual page with additional information on option
105153acd3b1SBarry Smith 
105253acd3b1SBarry Smith    Output Parameter:
105353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
105453acd3b1SBarry Smith 
105553acd3b1SBarry Smith    Level: intermediate
105653acd3b1SBarry Smith 
105753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
105853acd3b1SBarry Smith 
1059acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith     Concepts: options database^logical group
106253acd3b1SBarry Smith 
106353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1064acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
106553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
106653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1067acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
106853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
106953acd3b1SBarry Smith @*/
10707087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
107153acd3b1SBarry Smith {
107253acd3b1SBarry Smith   PetscErrorCode ierr;
10731ae3d29cSBarry Smith   PetscOptions   amsopt;
107453acd3b1SBarry Smith 
107553acd3b1SBarry Smith   PetscFunctionBegin;
10761ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10777781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1078ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1079a297a907SKarl Rupp 
1080ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10811ae3d29cSBarry Smith   }
108217326d04SJed Brown   *flg = PETSC_FALSE;
10830298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
108461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10852aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
108653acd3b1SBarry Smith   }
108753acd3b1SBarry Smith   PetscFunctionReturn(0);
108853acd3b1SBarry Smith }
108953acd3b1SBarry Smith 
109053acd3b1SBarry Smith #undef __FUNCT__
1091acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
109253acd3b1SBarry Smith /*@C
1093acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1094d5649816SBarry Smith        which at most a single value can be true.
109553acd3b1SBarry Smith 
10963f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
109753acd3b1SBarry Smith 
109853acd3b1SBarry Smith    Input Parameters:
109953acd3b1SBarry Smith +  opt - option name
110053acd3b1SBarry Smith .  text - short string that describes the option
110153acd3b1SBarry Smith -  man - manual page with additional information on option
110253acd3b1SBarry Smith 
110353acd3b1SBarry Smith    Output Parameter:
110453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
110553acd3b1SBarry Smith 
110653acd3b1SBarry Smith    Level: intermediate
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
110953acd3b1SBarry Smith 
1110acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
111153acd3b1SBarry Smith 
111253acd3b1SBarry Smith     Concepts: options database^logical group
111353acd3b1SBarry Smith 
111453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1115acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
111653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
111753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1118acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
111953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
112053acd3b1SBarry Smith @*/
11217087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
112253acd3b1SBarry Smith {
112353acd3b1SBarry Smith   PetscErrorCode ierr;
11241ae3d29cSBarry Smith   PetscOptions   amsopt;
112553acd3b1SBarry Smith 
112653acd3b1SBarry Smith   PetscFunctionBegin;
11271ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11287781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1129ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1130a297a907SKarl Rupp 
1131ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11321ae3d29cSBarry Smith   }
113317326d04SJed Brown   *flg = PETSC_FALSE;
11340298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
113561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11362aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
113753acd3b1SBarry Smith   }
113853acd3b1SBarry Smith   PetscFunctionReturn(0);
113953acd3b1SBarry Smith }
114053acd3b1SBarry Smith 
114153acd3b1SBarry Smith #undef __FUNCT__
1142acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
114353acd3b1SBarry Smith /*@C
1144acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
114553acd3b1SBarry Smith 
11463f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith    Input Parameters:
114953acd3b1SBarry Smith +  opt - option name
115053acd3b1SBarry Smith .  text - short string that describes the option
115153acd3b1SBarry Smith -  man - manual page with additional information on option
115253acd3b1SBarry Smith 
115353acd3b1SBarry Smith    Output Parameter:
115453acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
115553acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
115653acd3b1SBarry Smith 
115753acd3b1SBarry Smith    Level: beginner
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith    Concepts: options database^logical
116053acd3b1SBarry Smith 
116153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116253acd3b1SBarry Smith 
116353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1164acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1165acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
116653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
116753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1168acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
116953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
117053acd3b1SBarry Smith @*/
11717087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
117253acd3b1SBarry Smith {
117353acd3b1SBarry Smith   PetscErrorCode ierr;
1174ace3abfcSBarry Smith   PetscBool      iset;
1175aee2cecaSBarry Smith   PetscOptions   amsopt;
117653acd3b1SBarry Smith 
117753acd3b1SBarry Smith   PetscFunctionBegin;
1178b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
11797781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1180ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1181a297a907SKarl Rupp 
1182ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1183af6d86caSBarry Smith   }
1184acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
118553acd3b1SBarry Smith   if (!iset) {
118653acd3b1SBarry Smith     if (flg) *flg = deflt;
118753acd3b1SBarry Smith   }
118853acd3b1SBarry Smith   if (set) *set = iset;
118961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1190ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
11912aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
119253acd3b1SBarry Smith   }
119353acd3b1SBarry Smith   PetscFunctionReturn(0);
119453acd3b1SBarry Smith }
119553acd3b1SBarry Smith 
119653acd3b1SBarry Smith #undef __FUNCT__
119753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
119853acd3b1SBarry Smith /*@C
119953acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
120053acd3b1SBarry Smith    option in the database. The values must be separated with commas with
120153acd3b1SBarry Smith    no intervening spaces.
120253acd3b1SBarry Smith 
12033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120453acd3b1SBarry Smith 
120553acd3b1SBarry Smith    Input Parameters:
120653acd3b1SBarry Smith +  opt - the option one is seeking
120753acd3b1SBarry Smith .  text - short string describing option
120853acd3b1SBarry Smith .  man - manual page for option
120953acd3b1SBarry Smith -  nmax - maximum number of values
121053acd3b1SBarry Smith 
121153acd3b1SBarry Smith    Output Parameter:
121253acd3b1SBarry Smith +  value - location to copy values
121353acd3b1SBarry Smith .  nmax - actual number of values found
121453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
121553acd3b1SBarry Smith 
121653acd3b1SBarry Smith    Level: beginner
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith    Notes:
121953acd3b1SBarry Smith    The user should pass in an array of doubles
122053acd3b1SBarry Smith 
122153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122253acd3b1SBarry Smith 
122353acd3b1SBarry Smith    Concepts: options database^array of strings
122453acd3b1SBarry Smith 
122553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1226acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
122753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1229acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
123053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
123153acd3b1SBarry Smith @*/
12327087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
123353acd3b1SBarry Smith {
123453acd3b1SBarry Smith   PetscErrorCode ierr;
123553acd3b1SBarry Smith   PetscInt       i;
1236e26ddf31SBarry Smith   PetscOptions   amsopt;
123753acd3b1SBarry Smith 
123853acd3b1SBarry Smith   PetscFunctionBegin;
1239b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1240e26ddf31SBarry Smith     PetscReal *vals;
1241e26ddf31SBarry Smith 
1242e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1243e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1244e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1245e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1246e26ddf31SBarry Smith     amsopt->arraylength = *n;
1247e26ddf31SBarry Smith   }
124853acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
124961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1250a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
125153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1252a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
125353acd3b1SBarry Smith     }
12542aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
125553acd3b1SBarry Smith   }
125653acd3b1SBarry Smith   PetscFunctionReturn(0);
125753acd3b1SBarry Smith }
125853acd3b1SBarry Smith 
125953acd3b1SBarry Smith 
126053acd3b1SBarry Smith #undef __FUNCT__
126153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
126253acd3b1SBarry Smith /*@C
126353acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1264b32a342fSShri Abhyankar    option in the database.
126553acd3b1SBarry Smith 
12663f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126753acd3b1SBarry Smith 
126853acd3b1SBarry Smith    Input Parameters:
126953acd3b1SBarry Smith +  opt - the option one is seeking
127053acd3b1SBarry Smith .  text - short string describing option
127153acd3b1SBarry Smith .  man - manual page for option
1272f8a50e2bSBarry Smith -  n - maximum number of values
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Output Parameter:
127553acd3b1SBarry Smith +  value - location to copy values
1276f8a50e2bSBarry Smith .  n - actual number of values found
127753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Level: beginner
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith    Notes:
1282b32a342fSShri Abhyankar    The array can be passed as
1283b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12840fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12850fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12860fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1287b32a342fSShri Abhyankar 
1288b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
128953acd3b1SBarry Smith 
129053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129153acd3b1SBarry Smith 
1292b32a342fSShri Abhyankar    Concepts: options database^array of ints
129353acd3b1SBarry Smith 
129453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1295acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
129653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
129753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1298acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
129953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
130053acd3b1SBarry Smith @*/
13017087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
130253acd3b1SBarry Smith {
130353acd3b1SBarry Smith   PetscErrorCode ierr;
130453acd3b1SBarry Smith   PetscInt       i;
1305e26ddf31SBarry Smith   PetscOptions   amsopt;
130653acd3b1SBarry Smith 
130753acd3b1SBarry Smith   PetscFunctionBegin;
1308b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1309e26ddf31SBarry Smith     PetscInt *vals;
1310e26ddf31SBarry Smith 
1311e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1312e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1313e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1314e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1315e26ddf31SBarry Smith     amsopt->arraylength = *n;
1316e26ddf31SBarry Smith   }
131753acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
131861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
131953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
132053acd3b1SBarry Smith     for (i=1; i<*n; i++) {
132153acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
132253acd3b1SBarry Smith     }
13232aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
132453acd3b1SBarry Smith   }
132553acd3b1SBarry Smith   PetscFunctionReturn(0);
132653acd3b1SBarry Smith }
132753acd3b1SBarry Smith 
132853acd3b1SBarry Smith #undef __FUNCT__
132953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
133053acd3b1SBarry Smith /*@C
133153acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
133253acd3b1SBarry Smith    option in the database. The values must be separated with commas with
133353acd3b1SBarry Smith    no intervening spaces.
133453acd3b1SBarry Smith 
13353f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith    Input Parameters:
133853acd3b1SBarry Smith +  opt - the option one is seeking
133953acd3b1SBarry Smith .  text - short string describing option
134053acd3b1SBarry Smith .  man - manual page for option
134153acd3b1SBarry Smith -  nmax - maximum number of strings
134253acd3b1SBarry Smith 
134353acd3b1SBarry Smith    Output Parameter:
134453acd3b1SBarry Smith +  value - location to copy strings
134553acd3b1SBarry Smith .  nmax - actual number of strings found
134653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
134753acd3b1SBarry Smith 
134853acd3b1SBarry Smith    Level: beginner
134953acd3b1SBarry Smith 
135053acd3b1SBarry Smith    Notes:
135153acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
135253acd3b1SBarry Smith    strings returned by this function.
135353acd3b1SBarry Smith 
135453acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
135553acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
135853acd3b1SBarry Smith 
135953acd3b1SBarry Smith    Concepts: options database^array of strings
136053acd3b1SBarry Smith 
136153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1362acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
136353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
136453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1365acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
136653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
136753acd3b1SBarry Smith @*/
13687087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
136953acd3b1SBarry Smith {
137053acd3b1SBarry Smith   PetscErrorCode ierr;
13711ae3d29cSBarry Smith   PetscOptions   amsopt;
137253acd3b1SBarry Smith 
137353acd3b1SBarry Smith   PetscFunctionBegin;
13741ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13751ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13761ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1377a297a907SKarl Rupp 
13781ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13791ae3d29cSBarry Smith   }
138053acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
138161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13822aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
138353acd3b1SBarry Smith   }
138453acd3b1SBarry Smith   PetscFunctionReturn(0);
138553acd3b1SBarry Smith }
138653acd3b1SBarry Smith 
1387e2446a98SMatthew Knepley #undef __FUNCT__
1388acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1389e2446a98SMatthew Knepley /*@C
1390acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1391e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1392e2446a98SMatthew Knepley    no intervening spaces.
1393e2446a98SMatthew Knepley 
13943f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1395e2446a98SMatthew Knepley 
1396e2446a98SMatthew Knepley    Input Parameters:
1397e2446a98SMatthew Knepley +  opt - the option one is seeking
1398e2446a98SMatthew Knepley .  text - short string describing option
1399e2446a98SMatthew Knepley .  man - manual page for option
1400e2446a98SMatthew Knepley -  nmax - maximum number of values
1401e2446a98SMatthew Knepley 
1402e2446a98SMatthew Knepley    Output Parameter:
1403e2446a98SMatthew Knepley +  value - location to copy values
1404e2446a98SMatthew Knepley .  nmax - actual number of values found
1405e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1406e2446a98SMatthew Knepley 
1407e2446a98SMatthew Knepley    Level: beginner
1408e2446a98SMatthew Knepley 
1409e2446a98SMatthew Knepley    Notes:
1410e2446a98SMatthew Knepley    The user should pass in an array of doubles
1411e2446a98SMatthew Knepley 
1412e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1413e2446a98SMatthew Knepley 
1414e2446a98SMatthew Knepley    Concepts: options database^array of strings
1415e2446a98SMatthew Knepley 
1416e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1417acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1418e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1419e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1420acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1421e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1422e2446a98SMatthew Knepley @*/
14237087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1424e2446a98SMatthew Knepley {
1425e2446a98SMatthew Knepley   PetscErrorCode ierr;
1426e2446a98SMatthew Knepley   PetscInt       i;
14271ae3d29cSBarry Smith   PetscOptions   amsopt;
1428e2446a98SMatthew Knepley 
1429e2446a98SMatthew Knepley   PetscFunctionBegin;
14301ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1431ace3abfcSBarry Smith     PetscBool *vals;
14321ae3d29cSBarry Smith 
14337781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1434ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1435ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14361ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14371ae3d29cSBarry Smith     amsopt->arraylength = *n;
14381ae3d29cSBarry Smith   }
1439acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1440e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1441e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1442e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1443e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1444e2446a98SMatthew Knepley     }
14452aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1446e2446a98SMatthew Knepley   }
1447e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1448e2446a98SMatthew Knepley }
1449e2446a98SMatthew Knepley 
14508cc676e6SMatthew G Knepley #undef __FUNCT__
14518cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14528cc676e6SMatthew G Knepley /*@C
14538cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14548cc676e6SMatthew G Knepley 
14558cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14568cc676e6SMatthew G Knepley 
14578cc676e6SMatthew G Knepley    Input Parameters:
14588cc676e6SMatthew G Knepley +  opt - option name
14598cc676e6SMatthew G Knepley .  text - short string that describes the option
14608cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14618cc676e6SMatthew G Knepley 
14628cc676e6SMatthew G Knepley    Output Parameter:
14638cc676e6SMatthew G Knepley +  viewer - the viewer
14648cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14658cc676e6SMatthew G Knepley 
14668cc676e6SMatthew G Knepley    Level: beginner
14678cc676e6SMatthew G Knepley 
14688cc676e6SMatthew G Knepley    Concepts: options database^has int
14698cc676e6SMatthew G Knepley 
14708cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14718cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14728cc676e6SMatthew 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
14738cc676e6SMatthew G Knepley $                                     about the object to standard out
14748cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14758cc676e6SMatthew G Knepley $       draw
14768cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14778cc676e6SMatthew G Knepley 
1478cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14798cc676e6SMatthew G Knepley 
14808cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14818cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14828cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14838cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14848cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14858cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14868cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14878cc676e6SMatthew G Knepley @*/
1488cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14898cc676e6SMatthew G Knepley {
14908cc676e6SMatthew G Knepley   PetscErrorCode ierr;
14918cc676e6SMatthew G Knepley   PetscOptions   amsopt;
14928cc676e6SMatthew G Knepley 
14938cc676e6SMatthew G Knepley   PetscFunctionBegin;
14948cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
14958cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
14967781c08eSBarry Smith     amsopt->data = (void*)strdup("");
14978cc676e6SMatthew G Knepley   }
1498cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
14998cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15008cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15018cc676e6SMatthew G Knepley   }
15028cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15038cc676e6SMatthew G Knepley }
15048cc676e6SMatthew G Knepley 
150553acd3b1SBarry Smith 
150653acd3b1SBarry Smith #undef __FUNCT__
150753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
150853acd3b1SBarry Smith /*@C
1509b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
151053acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
151153acd3b1SBarry Smith 
15123f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
151353acd3b1SBarry Smith 
151453acd3b1SBarry Smith    Input Parameter:
151553acd3b1SBarry Smith .   head - the heading text
151653acd3b1SBarry Smith 
151753acd3b1SBarry Smith 
151853acd3b1SBarry Smith    Level: intermediate
151953acd3b1SBarry Smith 
152053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152153acd3b1SBarry Smith 
1522b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
152353acd3b1SBarry Smith 
152453acd3b1SBarry Smith    Concepts: options database^subheading
152553acd3b1SBarry Smith 
152653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1527acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
152853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
152953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1530acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
153153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
153253acd3b1SBarry Smith @*/
15337087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
153453acd3b1SBarry Smith {
153553acd3b1SBarry Smith   PetscErrorCode ierr;
153653acd3b1SBarry Smith 
153753acd3b1SBarry Smith   PetscFunctionBegin;
153861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
153953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
154053acd3b1SBarry Smith   }
154153acd3b1SBarry Smith   PetscFunctionReturn(0);
154253acd3b1SBarry Smith }
154353acd3b1SBarry Smith 
154453acd3b1SBarry Smith 
154553acd3b1SBarry Smith 
154653acd3b1SBarry Smith 
154753acd3b1SBarry Smith 
154853acd3b1SBarry Smith 
1549