xref: /petsc/src/sys/objects/aoptions.c (revision 64facd6c762c83bb0089ee53cdb0b0ed4101fc56)
17d0a6c19SBarry Smith 
253acd3b1SBarry Smith /*
33fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
43fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
553acd3b1SBarry Smith 
653acd3b1SBarry Smith */
753acd3b1SBarry Smith 
8afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
9665c2dedSJed Brown #include <petscviewer.h>
1053acd3b1SBarry Smith 
112aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
122aa6d131SJed Brown 
1353acd3b1SBarry Smith /*
1453acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
153fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1653acd3b1SBarry Smith 
1753acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1853acd3b1SBarry Smith */
19f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2053acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
2753acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
3253acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3353acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
346356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3753acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith 
410298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4253acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4361b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4553acd3b1SBarry Smith     }
4661b37b28SSatish Balay   }
4753acd3b1SBarry Smith   PetscFunctionReturn(0);
4853acd3b1SBarry Smith }
4953acd3b1SBarry Smith 
503194b578SJed Brown #undef __FUNCT__
513194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
553194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
563194b578SJed Brown {
573194b578SJed Brown   PetscErrorCode ierr;
583194b578SJed Brown   char           title[256];
593194b578SJed Brown   PetscBool      flg;
603194b578SJed Brown 
613194b578SJed Brown   PetscFunctionBegin;
623194b578SJed Brown   PetscValidHeader(obj,1);
633194b578SJed Brown   PetscOptionsObject.object         = obj;
643194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
65a297a907SKarl Rupp 
663194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
673194b578SJed Brown   if (flg) {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   } else {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   }
723194b578SJed Brown   ierr = PetscOptionsBegin_Private(obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
733194b578SJed Brown   PetscFunctionReturn(0);
743194b578SJed Brown }
753194b578SJed Brown 
7653acd3b1SBarry Smith /*
7753acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7853acd3b1SBarry Smith */
7953acd3b1SBarry Smith #undef __FUNCT__
8053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
81e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8253acd3b1SBarry Smith {
8353acd3b1SBarry Smith   int          ierr;
8453acd3b1SBarry Smith   PetscOptions next;
853be6e4c3SJed Brown   PetscBool    valid;
8653acd3b1SBarry Smith 
8753acd3b1SBarry Smith   PetscFunctionBegin;
883be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
893be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
903be6e4c3SJed Brown 
916e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9253acd3b1SBarry Smith   (*amsopt)->next = 0;
9353acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
946356e834SBarry Smith   (*amsopt)->type = t;
9553acd3b1SBarry Smith   (*amsopt)->data = 0;
9661b37b28SSatish Balay 
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9853acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
996356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10053acd3b1SBarry Smith 
101a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
102a297a907SKarl Rupp   else {
10353acd3b1SBarry Smith     next = PetscOptionsObject.next;
10453acd3b1SBarry Smith     while (next->next) next = next->next;
10553acd3b1SBarry Smith     next->next = *amsopt;
10653acd3b1SBarry Smith   }
10753acd3b1SBarry Smith   PetscFunctionReturn(0);
10853acd3b1SBarry Smith }
10953acd3b1SBarry Smith 
11053acd3b1SBarry Smith #undef __FUNCT__
111aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
112aee2cecaSBarry Smith /*
1133fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith     Collective on MPI_Comm
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith    Input Parameters:
1183fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1193fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1203fc1eb6aSBarry Smith -     str - location to store input
121aee2cecaSBarry Smith 
122aee2cecaSBarry Smith     Bugs:
123aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
124aee2cecaSBarry Smith 
125aee2cecaSBarry Smith */
1263fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
127aee2cecaSBarry Smith {
128330cf3c9SBarry Smith   size_t         i;
129aee2cecaSBarry Smith   char           c;
1303fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
131aee2cecaSBarry Smith   PetscErrorCode ierr;
132aee2cecaSBarry Smith 
133aee2cecaSBarry Smith   PetscFunctionBegin;
134aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
135aee2cecaSBarry Smith   if (!rank) {
136aee2cecaSBarry Smith     c = (char) getchar();
137aee2cecaSBarry Smith     i = 0;
138aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
139aee2cecaSBarry Smith       str[i++] = c;
140aee2cecaSBarry Smith       c = (char)getchar();
141aee2cecaSBarry Smith     }
142aee2cecaSBarry Smith     str[i] = 0;
143aee2cecaSBarry Smith   }
1444dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1453fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
146aee2cecaSBarry Smith   PetscFunctionReturn(0);
147aee2cecaSBarry Smith }
148aee2cecaSBarry Smith 
149aee2cecaSBarry Smith #undef __FUNCT__
150aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
151aee2cecaSBarry Smith /*
1523cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
153aee2cecaSBarry Smith 
154aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
155aee2cecaSBarry Smith 
1567781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1577781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1587781c08eSBarry Smith 
159aee2cecaSBarry Smith     Bugs:
1607781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1613cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
162aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
163aee2cecaSBarry Smith 
1643cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1653cc1e11dSBarry Smith      address space and communicating with the PETSc program
1663cc1e11dSBarry Smith 
167aee2cecaSBarry Smith */
168aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1696356e834SBarry Smith {
1706356e834SBarry Smith   PetscErrorCode ierr;
1716356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1726356e834SBarry Smith   char           str[512];
173a4404d99SBarry Smith   PetscInt       id;
1747781c08eSBarry Smith   PetscBool      bid;
175a4404d99SBarry Smith   PetscReal      ir,*valr;
176330cf3c9SBarry Smith   PetscInt       *vald;
177330cf3c9SBarry Smith   size_t         i;
1786356e834SBarry Smith 
179e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1806356e834SBarry Smith   while (next) {
1816356e834SBarry Smith     switch (next->type) {
1826356e834SBarry Smith     case OPTION_HEAD:
1836356e834SBarry Smith       break;
184e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
185e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
186e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
187e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
188e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
189e26ddf31SBarry Smith         if (i < next->arraylength-1) {
190e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
191e26ddf31SBarry Smith         }
192e26ddf31SBarry Smith       }
193e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
194e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
195e26ddf31SBarry Smith       if (str[0]) {
196e26ddf31SBarry Smith         PetscToken token;
197e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
198e26ddf31SBarry Smith         size_t     len;
199e26ddf31SBarry Smith         char       *value;
200ace3abfcSBarry Smith         PetscBool  foundrange;
201e26ddf31SBarry Smith 
202e26ddf31SBarry Smith         next->set = PETSC_TRUE;
203e26ddf31SBarry Smith         value     = str;
204e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
205e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
206e26ddf31SBarry Smith         while (n < nmax) {
207e26ddf31SBarry Smith           if (!value) break;
208e26ddf31SBarry Smith 
209e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
210e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
211e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
212e26ddf31SBarry Smith           if (value[0] == '-') i=2;
213e26ddf31SBarry Smith           else i=1;
214330cf3c9SBarry Smith           for (;i<len; i++) {
215e26ddf31SBarry Smith             if (value[i] == '-') {
216e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
217e26ddf31SBarry Smith               value[i] = 0;
218cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
219cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
220e32f2f54SBarry Smith               if (end <= start) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, %s-%s cannot have decreasing list",n,value,value+i+1);
221e32f2f54SBarry Smith               if (n + end - start - 1 >= nmax) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry, not enough space in left in array (%D) to contain entire range from %D to %D",n,nmax-n,start,end);
222e26ddf31SBarry Smith               for (; start<end; start++) {
223e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
224e26ddf31SBarry Smith               }
225e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
226e26ddf31SBarry Smith               break;
227e26ddf31SBarry Smith             }
228e26ddf31SBarry Smith           }
229e26ddf31SBarry Smith           if (!foundrange) {
230cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
231e26ddf31SBarry Smith             dvalue++;
232e26ddf31SBarry Smith             n++;
233e26ddf31SBarry Smith           }
234e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
235e26ddf31SBarry Smith         }
2368c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
237e26ddf31SBarry Smith       }
238e26ddf31SBarry Smith       break;
239e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
240e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
241e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
242e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
243e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
244e26ddf31SBarry Smith         if (i < next->arraylength-1) {
245e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
246e26ddf31SBarry Smith         }
247e26ddf31SBarry Smith       }
248e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
249e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
250e26ddf31SBarry Smith       if (str[0]) {
251e26ddf31SBarry Smith         PetscToken token;
252e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
253e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
254e26ddf31SBarry Smith         char       *value;
255e26ddf31SBarry Smith 
256e26ddf31SBarry Smith         next->set = PETSC_TRUE;
257e26ddf31SBarry Smith         value     = str;
258e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
259e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
260e26ddf31SBarry Smith         while (n < nmax) {
261e26ddf31SBarry Smith           if (!value) break;
262cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
263e26ddf31SBarry Smith           dvalue++;
264e26ddf31SBarry Smith           n++;
265e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
266e26ddf31SBarry Smith         }
2678c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
268e26ddf31SBarry Smith       }
269e26ddf31SBarry Smith       break;
2706356e834SBarry Smith     case OPTION_INT:
271e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%d>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(int*)next->data,next->text,next->man);CHKERRQ(ierr);
2723fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2733fc1eb6aSBarry Smith       if (str[0]) {
274c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
275c272547aSJed Brown         sscanf(str,"%lld",&id);
276c272547aSJed Brown #else
277aee2cecaSBarry Smith         sscanf(str,"%d",&id);
278c272547aSJed Brown #endif
279aee2cecaSBarry Smith         next->set = PETSC_TRUE;
280a297a907SKarl Rupp 
281aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
282aee2cecaSBarry Smith       }
283aee2cecaSBarry Smith       break;
284aee2cecaSBarry Smith     case OPTION_REAL:
285e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%g>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(double*)next->data,next->text,next->man);CHKERRQ(ierr);
2863fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2873fc1eb6aSBarry Smith       if (str[0]) {
288ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
289a4404d99SBarry Smith         sscanf(str,"%e",&ir);
290ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
291aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
292ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
293d9822059SBarry Smith         ir = strtoflt128(str,0);
294d9822059SBarry Smith #else
295513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
296a4404d99SBarry Smith #endif
297aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
298aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
299aee2cecaSBarry Smith       }
300aee2cecaSBarry Smith       break;
3017781c08eSBarry Smith     case OPTION_BOOL:
3027781c08eSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3037781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3047781c08eSBarry Smith       if (str[0]) {
3057781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3067781c08eSBarry Smith         next->set = PETSC_TRUE;
3077781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3087781c08eSBarry Smith       }
3097781c08eSBarry Smith       break;
310aee2cecaSBarry Smith     case OPTION_STRING:
311e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,(char*)next->data,next->text,next->man);CHKERRQ(ierr);
3123fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3133fc1eb6aSBarry Smith       if (str[0]) {
314aee2cecaSBarry Smith         next->set = PETSC_TRUE;
315*64facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
316c979a496SBarry Smith         next->data = (void*)strdup(str);
3176356e834SBarry Smith       }
3186356e834SBarry Smith       break;
319a264d7a6SBarry Smith     case OPTION_FLIST:
320140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3213cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3223cc1e11dSBarry Smith       if (str[0]) {
3233cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3243cc1e11dSBarry Smith         next->set = PETSC_TRUE;
325*64facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
326c979a496SBarry Smith         next->data = (void*)strdup(str);
3273cc1e11dSBarry Smith       }
3283cc1e11dSBarry Smith       break;
329b432afdaSMatthew Knepley     default:
330b432afdaSMatthew Knepley       break;
3316356e834SBarry Smith     }
3326356e834SBarry Smith     next = next->next;
3336356e834SBarry Smith   }
3346356e834SBarry Smith   PetscFunctionReturn(0);
3356356e834SBarry Smith }
3366356e834SBarry Smith 
337e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
338e04113cfSBarry Smith #include <petscviewersaws.h>
339d5649816SBarry Smith 
340d5649816SBarry Smith static int count = 0;
341d5649816SBarry Smith 
342b3506946SBarry Smith #undef __FUNCT__
343e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
344e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
345d5649816SBarry Smith {
3462657e9d9SBarry Smith   PetscFunctionBegin;
347d5649816SBarry Smith   PetscFunctionReturn(0);
348d5649816SBarry Smith }
349d5649816SBarry Smith 
350d5649816SBarry Smith #undef __FUNCT__
3517781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
352b3506946SBarry Smith /*
3537781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
354b3506946SBarry Smith 
355b3506946SBarry Smith     Bugs:
356b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
357b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
358b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
359b3506946SBarry Smith 
360b3506946SBarry Smith 
361b3506946SBarry Smith */
362475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
363b3506946SBarry Smith {
364b3506946SBarry Smith   PetscErrorCode ierr;
365b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
366d5649816SBarry Smith   static int     mancount = 0;
367b3506946SBarry Smith   char           options[16];
3687aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
36988a9752cSBarry Smith   char           manname[16],textname[16];
3702657e9d9SBarry Smith   char           dir[1024];
371b3506946SBarry Smith 
372b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
373b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
374a297a907SKarl Rupp 
375e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
3761bc75a8dSBarry Smith 
377d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
3787781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
3797781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
3802657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
3812657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
382b3506946SBarry Smith 
383b3506946SBarry Smith   while (next) {
384d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
3852657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
3867781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
387d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
3887781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
3897781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
3909f32e415SBarry Smith 
391b3506946SBarry Smith     switch (next->type) {
392b3506946SBarry Smith     case OPTION_HEAD:
393b3506946SBarry Smith       break;
394b3506946SBarry Smith     case OPTION_INT_ARRAY:
3957781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
3962657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
397b3506946SBarry Smith       break;
398b3506946SBarry Smith     case OPTION_REAL_ARRAY:
3997781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4002657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
401b3506946SBarry Smith       break;
402b3506946SBarry Smith     case OPTION_INT:
4037781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4042657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
405b3506946SBarry Smith       break;
406b3506946SBarry Smith     case OPTION_REAL:
4077781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4082657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
409b3506946SBarry Smith       break;
4107781c08eSBarry Smith     case OPTION_BOOL:
4117781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4122657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4131ae3d29cSBarry Smith       break;
4147781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4157781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4162657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
41771f08665SBarry Smith       break;
418b3506946SBarry Smith     case OPTION_STRING:
4197781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4207781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4211ae3d29cSBarry Smith       break;
4221ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4237781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4242657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
425b3506946SBarry Smith       break;
426a264d7a6SBarry Smith     case OPTION_FLIST:
427a264d7a6SBarry Smith       {
428a264d7a6SBarry 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));
431a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
432a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
433a264d7a6SBarry Smith       }
434a264d7a6SBarry Smith       break;
4351ae3d29cSBarry Smith     case OPTION_ELIST:
436a264d7a6SBarry Smith       {
437a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4387781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4397781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
440a264d7a6SBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char*),&next->edata);CHKERRQ(ierr);
4411ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
442a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
443a264d7a6SBarry Smith       }
444a264d7a6SBarry Smith       break;
445b3506946SBarry Smith     default:
446b3506946SBarry Smith       break;
447b3506946SBarry Smith     }
448b3506946SBarry Smith     next = next->next;
449b3506946SBarry Smith   }
450b3506946SBarry Smith 
451b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4527aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
453b3506946SBarry Smith 
45488a9752cSBarry Smith   /* determine if any values have been set in GUI */
45588a9752cSBarry Smith   next = PetscOptionsObject.next;
45688a9752cSBarry Smith   while (next) {
45788a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
45888a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
45988a9752cSBarry Smith     next = next->next;
46088a9752cSBarry Smith   }
46188a9752cSBarry Smith 
462b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
463b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
464b3506946SBarry Smith 
4659a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
466b3506946SBarry Smith   PetscFunctionReturn(0);
467b3506946SBarry Smith }
468b3506946SBarry Smith #endif
469b3506946SBarry Smith 
4706356e834SBarry Smith #undef __FUNCT__
47153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
47253acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47353acd3b1SBarry Smith {
47453acd3b1SBarry Smith   PetscErrorCode ierr;
4756356e834SBarry Smith   PetscOptions   last;
4766356e834SBarry Smith   char           option[256],value[1024],tmp[32];
477330cf3c9SBarry Smith   size_t         j;
47853acd3b1SBarry Smith 
47953acd3b1SBarry Smith   PetscFunctionBegin;
480aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
481b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
482a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
483475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
484b3506946SBarry Smith #else
48571f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
486b3506946SBarry Smith #endif
487aee2cecaSBarry Smith     }
488aee2cecaSBarry Smith   }
4896356e834SBarry Smith 
490c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
491c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4926356e834SBarry Smith 
4936356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4946356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49561b37b28SSatish Balay   /* reset alreadyprinted flag */
49661b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4973194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4980298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4996356e834SBarry Smith 
5006356e834SBarry Smith   while (PetscOptionsObject.next) {
5016356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5026356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5036356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5046356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5056356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5066356e834SBarry Smith       } else {
5076356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5086356e834SBarry Smith       }
5096356e834SBarry Smith 
5106356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5116356e834SBarry Smith       case OPTION_HEAD:
5126356e834SBarry Smith         break;
513e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
514e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
515e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
516e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
517e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
518e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
519e26ddf31SBarry Smith         }
520e26ddf31SBarry Smith         break;
5216356e834SBarry Smith       case OPTION_INT:
5227a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5236356e834SBarry Smith         break;
5246356e834SBarry Smith       case OPTION_REAL:
5257a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5266356e834SBarry Smith         break;
5276356e834SBarry Smith       case OPTION_REAL_ARRAY:
5287a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5296356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5307a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5316356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5326356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5336356e834SBarry Smith         }
5346356e834SBarry Smith         break;
5357781c08eSBarry Smith       case OPTION_BOOL:
53671f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5376356e834SBarry Smith         break;
5387781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
539ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5401ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
541ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5421ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5431ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5441ae3d29cSBarry Smith         }
5451ae3d29cSBarry Smith         break;
546a264d7a6SBarry Smith       case OPTION_FLIST:
5471ae3d29cSBarry Smith       case OPTION_ELIST:
5487781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5496356e834SBarry Smith         break;
5501ae3d29cSBarry Smith       case OPTION_STRING:
551475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5528da2146eSBarry Smith         break;
5531ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5541ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5551ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5561ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5571ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5581ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5591ae3d29cSBarry Smith         }
5606356e834SBarry Smith         break;
5616356e834SBarry Smith       }
5626356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5636356e834SBarry Smith     }
564503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
565503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5666356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56771f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
568c979a496SBarry Smith 
569c979a496SBarry Smith     if ((PetscOptionsObject.next->type == OPTION_STRING) || (PetscOptionsObject.next->type == OPTION_FLIST) || (PetscOptionsObject.next->type == OPTION_ELIST)){
570*64facd6cSBarry Smith       /* must use system free since SAWs may have allocated it */
571c979a496SBarry Smith       free(PetscOptionsObject.next->data);
572c979a496SBarry Smith     } else {
5737781c08eSBarry Smith       ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
574c979a496SBarry Smith     }
5757781c08eSBarry Smith 
5766356e834SBarry Smith     last                    = PetscOptionsObject.next;
5776356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5786356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5796356e834SBarry Smith   }
5806356e834SBarry Smith   PetscOptionsObject.next = 0;
58153acd3b1SBarry Smith   PetscFunctionReturn(0);
58253acd3b1SBarry Smith }
58353acd3b1SBarry Smith 
58453acd3b1SBarry Smith #undef __FUNCT__
58553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
58653acd3b1SBarry Smith /*@C
58753acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58853acd3b1SBarry Smith 
5893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
59053acd3b1SBarry Smith 
59153acd3b1SBarry Smith    Input Parameters:
59253acd3b1SBarry Smith +  opt - option name
59353acd3b1SBarry Smith .  text - short string that describes the option
59453acd3b1SBarry Smith .  man - manual page with additional information on option
59553acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
59653acd3b1SBarry Smith -  defaultv - the default (current) value
59753acd3b1SBarry Smith 
59853acd3b1SBarry Smith    Output Parameter:
59953acd3b1SBarry Smith +  value - the  value to return
600b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
60153acd3b1SBarry Smith 
60253acd3b1SBarry Smith    Level: beginner
60353acd3b1SBarry Smith 
60453acd3b1SBarry Smith    Concepts: options database
60553acd3b1SBarry Smith 
60653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
60753acd3b1SBarry Smith 
60853acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60953acd3b1SBarry Smith 
61053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
611acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
612acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
61353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
61453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
615acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
616a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
61753acd3b1SBarry Smith @*/
6187087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61953acd3b1SBarry Smith {
62053acd3b1SBarry Smith   PetscErrorCode ierr;
62153acd3b1SBarry Smith   PetscInt       ntext = 0;
622aa5bb8c0SSatish Balay   PetscInt       tval;
623ace3abfcSBarry Smith   PetscBool      tflg;
62453acd3b1SBarry Smith 
62553acd3b1SBarry Smith   PetscFunctionBegin;
62653acd3b1SBarry Smith   while (list[ntext++]) {
627e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62853acd3b1SBarry Smith   }
629e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
63053acd3b1SBarry Smith   ntext -= 3;
631aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
632aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
633aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
634aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
63553acd3b1SBarry Smith   PetscFunctionReturn(0);
63653acd3b1SBarry Smith }
63753acd3b1SBarry Smith 
63853acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63953acd3b1SBarry Smith #undef __FUNCT__
64053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
64153acd3b1SBarry Smith /*@C
64253acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
64353acd3b1SBarry Smith 
6443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Input Parameters:
64753acd3b1SBarry Smith +  opt - option name
64853acd3b1SBarry Smith .  text - short string that describes the option
64953acd3b1SBarry Smith .  man - manual page with additional information on option
65053acd3b1SBarry Smith -  defaultv - the default (current) value
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith    Output Parameter:
65353acd3b1SBarry Smith +  value - the integer value to return
65453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
65553acd3b1SBarry Smith 
65653acd3b1SBarry Smith    Level: beginner
65753acd3b1SBarry Smith 
65853acd3b1SBarry Smith    Concepts: options database^has int
65953acd3b1SBarry Smith 
66053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
66153acd3b1SBarry Smith 
66253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
663acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
664acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
667acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
668a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66953acd3b1SBarry Smith @*/
6707087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
67153acd3b1SBarry Smith {
67253acd3b1SBarry Smith   PetscErrorCode ierr;
6736356e834SBarry Smith   PetscOptions   amsopt;
67453acd3b1SBarry Smith 
67553acd3b1SBarry Smith   PetscFunctionBegin;
676b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6776356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6786356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
679a297a907SKarl Rupp 
6806356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
681af6d86caSBarry Smith   }
68253acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
68361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6842aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
68553acd3b1SBarry Smith   }
68653acd3b1SBarry Smith   PetscFunctionReturn(0);
68753acd3b1SBarry Smith }
68853acd3b1SBarry Smith 
68953acd3b1SBarry Smith #undef __FUNCT__
69053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
69153acd3b1SBarry Smith /*@C
69253acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
69353acd3b1SBarry Smith 
6943f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith    Input Parameters:
69753acd3b1SBarry Smith +  opt - option name
69853acd3b1SBarry Smith .  text - short string that describes the option
69953acd3b1SBarry Smith .  man - manual page with additional information on option
700bcbf2dc5SJed Brown .  defaultv - the default (current) value
701bcbf2dc5SJed Brown -  len - length of the result string including null terminator
70253acd3b1SBarry Smith 
70353acd3b1SBarry Smith    Output Parameter:
70453acd3b1SBarry Smith +  value - the value to return
70553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70653acd3b1SBarry Smith 
70753acd3b1SBarry Smith    Level: beginner
70853acd3b1SBarry Smith 
70953acd3b1SBarry Smith    Concepts: options database^has int
71053acd3b1SBarry Smith 
71153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
71253acd3b1SBarry Smith 
7137fccdfe4SBarry 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).
7147fccdfe4SBarry Smith 
71553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
716acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
717acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
720acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
721a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
72253acd3b1SBarry Smith @*/
7237087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
72453acd3b1SBarry Smith {
72553acd3b1SBarry Smith   PetscErrorCode ierr;
726aee2cecaSBarry Smith   PetscOptions   amsopt;
72753acd3b1SBarry Smith 
72853acd3b1SBarry Smith   PetscFunctionBegin;
729b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
730aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
731*64facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
732c979a496SBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
733af6d86caSBarry Smith   }
73453acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
73561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7362aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73753acd3b1SBarry Smith   }
73853acd3b1SBarry Smith   PetscFunctionReturn(0);
73953acd3b1SBarry Smith }
74053acd3b1SBarry Smith 
74153acd3b1SBarry Smith #undef __FUNCT__
74253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
74353acd3b1SBarry Smith /*@C
74453acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
74553acd3b1SBarry Smith 
7463f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74753acd3b1SBarry Smith 
74853acd3b1SBarry Smith    Input Parameters:
74953acd3b1SBarry Smith +  opt - option name
75053acd3b1SBarry Smith .  text - short string that describes the option
75153acd3b1SBarry Smith .  man - manual page with additional information on option
75253acd3b1SBarry Smith -  defaultv - the default (current) value
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith    Output Parameter:
75553acd3b1SBarry Smith +  value - the value to return
75653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith    Level: beginner
75953acd3b1SBarry Smith 
76053acd3b1SBarry Smith    Concepts: options database^has int
76153acd3b1SBarry Smith 
76253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
76353acd3b1SBarry Smith 
76453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
765acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
766acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
76753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
769acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
770a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
77153acd3b1SBarry Smith @*/
7727087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
77353acd3b1SBarry Smith {
77453acd3b1SBarry Smith   PetscErrorCode ierr;
775538aa990SBarry Smith   PetscOptions   amsopt;
77653acd3b1SBarry Smith 
77753acd3b1SBarry Smith   PetscFunctionBegin;
778b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
779538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
780538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
781a297a907SKarl Rupp 
782538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
783538aa990SBarry Smith   }
78453acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
78561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7862aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78753acd3b1SBarry Smith   }
78853acd3b1SBarry Smith   PetscFunctionReturn(0);
78953acd3b1SBarry Smith }
79053acd3b1SBarry Smith 
79153acd3b1SBarry Smith #undef __FUNCT__
79253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
79353acd3b1SBarry Smith /*@C
79453acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
79553acd3b1SBarry Smith 
7963f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith    Input Parameters:
79953acd3b1SBarry Smith +  opt - option name
80053acd3b1SBarry Smith .  text - short string that describes the option
80153acd3b1SBarry Smith .  man - manual page with additional information on option
80253acd3b1SBarry Smith -  defaultv - the default (current) value
80353acd3b1SBarry Smith 
80453acd3b1SBarry Smith    Output Parameter:
80553acd3b1SBarry Smith +  value - the value to return
80653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80753acd3b1SBarry Smith 
80853acd3b1SBarry Smith    Level: beginner
80953acd3b1SBarry Smith 
81053acd3b1SBarry Smith    Concepts: options database^has int
81153acd3b1SBarry Smith 
81253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
81353acd3b1SBarry Smith 
81453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
815acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
816acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
819acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
820a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
82153acd3b1SBarry Smith @*/
8227087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
82353acd3b1SBarry Smith {
82453acd3b1SBarry Smith   PetscErrorCode ierr;
82553acd3b1SBarry Smith 
82653acd3b1SBarry Smith   PetscFunctionBegin;
82753acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
82853acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82953acd3b1SBarry Smith #else
83053acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
83153acd3b1SBarry Smith #endif
83253acd3b1SBarry Smith   PetscFunctionReturn(0);
83353acd3b1SBarry Smith }
83453acd3b1SBarry Smith 
83553acd3b1SBarry Smith #undef __FUNCT__
83653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
83753acd3b1SBarry Smith /*@C
83890d69ab7SBarry 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
83990d69ab7SBarry Smith                       its value is set to false.
84053acd3b1SBarry Smith 
8413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith    Input Parameters:
84453acd3b1SBarry Smith +  opt - option name
84553acd3b1SBarry Smith .  text - short string that describes the option
84653acd3b1SBarry Smith -  man - manual page with additional information on option
84753acd3b1SBarry Smith 
84853acd3b1SBarry Smith    Output Parameter:
84953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Level: beginner
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith    Concepts: options database^has int
85453acd3b1SBarry Smith 
85553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
858acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
859acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
86053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
86153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
862acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
863a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86453acd3b1SBarry Smith @*/
8657087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
86653acd3b1SBarry Smith {
86753acd3b1SBarry Smith   PetscErrorCode ierr;
8681ae3d29cSBarry Smith   PetscOptions   amsopt;
86953acd3b1SBarry Smith 
87053acd3b1SBarry Smith   PetscFunctionBegin;
8711ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8727781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
873ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
874a297a907SKarl Rupp 
875ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8761ae3d29cSBarry Smith   }
87753acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
87861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8792aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
88053acd3b1SBarry Smith   }
88153acd3b1SBarry Smith   PetscFunctionReturn(0);
88253acd3b1SBarry Smith }
88353acd3b1SBarry Smith 
88453acd3b1SBarry Smith #undef __FUNCT__
885a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
88653acd3b1SBarry Smith /*@C
887a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
88853acd3b1SBarry Smith 
8893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
89053acd3b1SBarry Smith 
89153acd3b1SBarry Smith    Input Parameters:
89253acd3b1SBarry Smith +  opt - option name
89353acd3b1SBarry Smith .  text - short string that describes the option
89453acd3b1SBarry Smith .  man - manual page with additional information on option
89553acd3b1SBarry Smith .  list - the possible choices
8963cc1e11dSBarry Smith .  defaultv - the default (current) value
8973cc1e11dSBarry Smith -  len - the length of the character array value
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith    Output Parameter:
90053acd3b1SBarry Smith +  value - the value to return
90153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
90253acd3b1SBarry Smith 
90353acd3b1SBarry Smith    Level: intermediate
90453acd3b1SBarry Smith 
90553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
90653acd3b1SBarry Smith 
90753acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
90853acd3b1SBarry Smith 
90953acd3b1SBarry Smith    To get a listing of all currently specified options,
91088c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
91153acd3b1SBarry Smith 
912eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
913eabe10d7SBarry Smith 
91453acd3b1SBarry Smith    Concepts: options database^list
91553acd3b1SBarry Smith 
91653acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
917acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
91853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
920acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
921a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
92253acd3b1SBarry Smith @*/
923a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
92453acd3b1SBarry Smith {
92553acd3b1SBarry Smith   PetscErrorCode ierr;
9263cc1e11dSBarry Smith   PetscOptions   amsopt;
92753acd3b1SBarry Smith 
92853acd3b1SBarry Smith   PetscFunctionBegin;
929b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
930a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
931*64facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
932c979a496SBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9333cc1e11dSBarry Smith     amsopt->flist = list;
9343cc1e11dSBarry Smith   }
93553acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
93661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
937140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
93853acd3b1SBarry Smith   }
93953acd3b1SBarry Smith   PetscFunctionReturn(0);
94053acd3b1SBarry Smith }
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith #undef __FUNCT__
94353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
94453acd3b1SBarry Smith /*@C
94553acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
94653acd3b1SBarry Smith 
9473f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94853acd3b1SBarry Smith 
94953acd3b1SBarry Smith    Input Parameters:
95053acd3b1SBarry Smith +  opt - option name
95153acd3b1SBarry Smith .  ltext - short string that describes the option
95253acd3b1SBarry Smith .  man - manual page with additional information on option
953a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
95453acd3b1SBarry Smith .  ntext - number of choices
95553acd3b1SBarry Smith -  defaultv - the default (current) value
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith    Output Parameter:
95853acd3b1SBarry Smith +  value - the index of the value to return
95953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith    Level: intermediate
96253acd3b1SBarry Smith 
96353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
96453acd3b1SBarry Smith 
965a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
96653acd3b1SBarry Smith 
96753acd3b1SBarry Smith    Concepts: options database^list
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
970acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
97153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
97253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
973acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
974a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
97553acd3b1SBarry Smith @*/
9767087cfbeSBarry 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)
97753acd3b1SBarry Smith {
97853acd3b1SBarry Smith   PetscErrorCode ierr;
97953acd3b1SBarry Smith   PetscInt       i;
9801ae3d29cSBarry Smith   PetscOptions   amsopt;
98153acd3b1SBarry Smith 
98253acd3b1SBarry Smith   PetscFunctionBegin;
9831ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
984d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
985*64facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
986c979a496SBarry Smith     amsopt->data = (void*)strdup(defaultv ? defaultv : "");
9871ae3d29cSBarry Smith     amsopt->list  = list;
9881ae3d29cSBarry Smith     amsopt->nlist = ntext;
9891ae3d29cSBarry Smith   }
99053acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
99161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
99253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
99353acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
99453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
99553acd3b1SBarry Smith     }
9962aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
99753acd3b1SBarry Smith   }
99853acd3b1SBarry Smith   PetscFunctionReturn(0);
99953acd3b1SBarry Smith }
100053acd3b1SBarry Smith 
100153acd3b1SBarry Smith #undef __FUNCT__
1002acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
100353acd3b1SBarry Smith /*@C
1004acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1005d5649816SBarry Smith        which at most a single value can be true.
100653acd3b1SBarry Smith 
10073f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Input Parameters:
101053acd3b1SBarry Smith +  opt - option name
101153acd3b1SBarry Smith .  text - short string that describes the option
101253acd3b1SBarry Smith -  man - manual page with additional information on option
101353acd3b1SBarry Smith 
101453acd3b1SBarry Smith    Output Parameter:
101553acd3b1SBarry Smith .  flg - whether that option was set or not
101653acd3b1SBarry Smith 
101753acd3b1SBarry Smith    Level: intermediate
101853acd3b1SBarry Smith 
101953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
102053acd3b1SBarry Smith 
1021acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
102253acd3b1SBarry Smith 
102353acd3b1SBarry Smith     Concepts: options database^logical group
102453acd3b1SBarry Smith 
102553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1026acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1029acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1030a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
103153acd3b1SBarry Smith @*/
10327087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
103353acd3b1SBarry Smith {
103453acd3b1SBarry Smith   PetscErrorCode ierr;
10351ae3d29cSBarry Smith   PetscOptions   amsopt;
103653acd3b1SBarry Smith 
103753acd3b1SBarry Smith   PetscFunctionBegin;
10381ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10397781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1040ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1041a297a907SKarl Rupp 
1042ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10431ae3d29cSBarry Smith   }
104468b16fdaSBarry Smith   *flg = PETSC_FALSE;
10450298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
104661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10482aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104953acd3b1SBarry Smith   }
105053acd3b1SBarry Smith   PetscFunctionReturn(0);
105153acd3b1SBarry Smith }
105253acd3b1SBarry Smith 
105353acd3b1SBarry Smith #undef __FUNCT__
1054acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
105553acd3b1SBarry Smith /*@C
1056acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1057d5649816SBarry Smith        which at most a single value can be true.
105853acd3b1SBarry Smith 
10593f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Input Parameters:
106253acd3b1SBarry Smith +  opt - option name
106353acd3b1SBarry Smith .  text - short string that describes the option
106453acd3b1SBarry Smith -  man - manual page with additional information on option
106553acd3b1SBarry Smith 
106653acd3b1SBarry Smith    Output Parameter:
106753acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106853acd3b1SBarry Smith 
106953acd3b1SBarry Smith    Level: intermediate
107053acd3b1SBarry Smith 
107153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
107253acd3b1SBarry Smith 
1073acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
107453acd3b1SBarry Smith 
107553acd3b1SBarry Smith     Concepts: options database^logical group
107653acd3b1SBarry Smith 
107753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1078acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
108053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1081acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1082a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
108353acd3b1SBarry Smith @*/
10847087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
108553acd3b1SBarry Smith {
108653acd3b1SBarry Smith   PetscErrorCode ierr;
10871ae3d29cSBarry Smith   PetscOptions   amsopt;
108853acd3b1SBarry Smith 
108953acd3b1SBarry Smith   PetscFunctionBegin;
10901ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10917781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1092ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1093a297a907SKarl Rupp 
1094ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10951ae3d29cSBarry Smith   }
109617326d04SJed Brown   *flg = PETSC_FALSE;
10970298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10992aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
110053acd3b1SBarry Smith   }
110153acd3b1SBarry Smith   PetscFunctionReturn(0);
110253acd3b1SBarry Smith }
110353acd3b1SBarry Smith 
110453acd3b1SBarry Smith #undef __FUNCT__
1105acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
110653acd3b1SBarry Smith /*@C
1107acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1108d5649816SBarry Smith        which at most a single value can be true.
110953acd3b1SBarry Smith 
11103f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
111153acd3b1SBarry Smith 
111253acd3b1SBarry Smith    Input Parameters:
111353acd3b1SBarry Smith +  opt - option name
111453acd3b1SBarry Smith .  text - short string that describes the option
111553acd3b1SBarry Smith -  man - manual page with additional information on option
111653acd3b1SBarry Smith 
111753acd3b1SBarry Smith    Output Parameter:
111853acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111953acd3b1SBarry Smith 
112053acd3b1SBarry Smith    Level: intermediate
112153acd3b1SBarry Smith 
112253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
112353acd3b1SBarry Smith 
1124acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
112553acd3b1SBarry Smith 
112653acd3b1SBarry Smith     Concepts: options database^logical group
112753acd3b1SBarry Smith 
112853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1129acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
113053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
113153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1132acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1133a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
113453acd3b1SBarry Smith @*/
11357087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
113653acd3b1SBarry Smith {
113753acd3b1SBarry Smith   PetscErrorCode ierr;
11381ae3d29cSBarry Smith   PetscOptions   amsopt;
113953acd3b1SBarry Smith 
114053acd3b1SBarry Smith   PetscFunctionBegin;
11411ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11427781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1143ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1144a297a907SKarl Rupp 
1145ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11461ae3d29cSBarry Smith   }
114717326d04SJed Brown   *flg = PETSC_FALSE;
11480298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11502aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
115153acd3b1SBarry Smith   }
115253acd3b1SBarry Smith   PetscFunctionReturn(0);
115353acd3b1SBarry Smith }
115453acd3b1SBarry Smith 
115553acd3b1SBarry Smith #undef __FUNCT__
1156acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
115753acd3b1SBarry Smith /*@C
1158acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115953acd3b1SBarry Smith 
11603f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
116153acd3b1SBarry Smith 
116253acd3b1SBarry Smith    Input Parameters:
116353acd3b1SBarry Smith +  opt - option name
116453acd3b1SBarry Smith .  text - short string that describes the option
116553acd3b1SBarry Smith -  man - manual page with additional information on option
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith    Output Parameter:
116853acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116953acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
117053acd3b1SBarry Smith 
117153acd3b1SBarry Smith    Level: beginner
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith    Concepts: options database^logical
117453acd3b1SBarry Smith 
117553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117653acd3b1SBarry Smith 
117753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1178acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1179acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
118053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
118153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1182acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1183a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
118453acd3b1SBarry Smith @*/
11857087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
118653acd3b1SBarry Smith {
118753acd3b1SBarry Smith   PetscErrorCode ierr;
1188ace3abfcSBarry Smith   PetscBool      iset;
1189aee2cecaSBarry Smith   PetscOptions   amsopt;
119053acd3b1SBarry Smith 
119153acd3b1SBarry Smith   PetscFunctionBegin;
1192b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
11937781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1194ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1195a297a907SKarl Rupp 
1196ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1197af6d86caSBarry Smith   }
1198acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119953acd3b1SBarry Smith   if (!iset) {
120053acd3b1SBarry Smith     if (flg) *flg = deflt;
120153acd3b1SBarry Smith   }
120253acd3b1SBarry Smith   if (set) *set = iset;
120361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1204ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12052aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
120653acd3b1SBarry Smith   }
120753acd3b1SBarry Smith   PetscFunctionReturn(0);
120853acd3b1SBarry Smith }
120953acd3b1SBarry Smith 
121053acd3b1SBarry Smith #undef __FUNCT__
121153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
121253acd3b1SBarry Smith /*@C
121353acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
121453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
121553acd3b1SBarry Smith    no intervening spaces.
121653acd3b1SBarry Smith 
12173f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121853acd3b1SBarry Smith 
121953acd3b1SBarry Smith    Input Parameters:
122053acd3b1SBarry Smith +  opt - the option one is seeking
122153acd3b1SBarry Smith .  text - short string describing option
122253acd3b1SBarry Smith .  man - manual page for option
122353acd3b1SBarry Smith -  nmax - maximum number of values
122453acd3b1SBarry Smith 
122553acd3b1SBarry Smith    Output Parameter:
122653acd3b1SBarry Smith +  value - location to copy values
122753acd3b1SBarry Smith .  nmax - actual number of values found
122853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122953acd3b1SBarry Smith 
123053acd3b1SBarry Smith    Level: beginner
123153acd3b1SBarry Smith 
123253acd3b1SBarry Smith    Notes:
123353acd3b1SBarry Smith    The user should pass in an array of doubles
123453acd3b1SBarry Smith 
123553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith    Concepts: options database^array of strings
123853acd3b1SBarry Smith 
123953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1240acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
124153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
124253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1243acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1244a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
124553acd3b1SBarry Smith @*/
12467087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
124753acd3b1SBarry Smith {
124853acd3b1SBarry Smith   PetscErrorCode ierr;
124953acd3b1SBarry Smith   PetscInt       i;
1250e26ddf31SBarry Smith   PetscOptions   amsopt;
125153acd3b1SBarry Smith 
125253acd3b1SBarry Smith   PetscFunctionBegin;
1253b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1254e26ddf31SBarry Smith     PetscReal *vals;
1255e26ddf31SBarry Smith 
1256e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1257e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1258e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1259e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1260e26ddf31SBarry Smith     amsopt->arraylength = *n;
1261e26ddf31SBarry Smith   }
126253acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
126361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1264a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
126553acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1266a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
126753acd3b1SBarry Smith     }
12682aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126953acd3b1SBarry Smith   }
127053acd3b1SBarry Smith   PetscFunctionReturn(0);
127153acd3b1SBarry Smith }
127253acd3b1SBarry Smith 
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith #undef __FUNCT__
127553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
127653acd3b1SBarry Smith /*@C
127753acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1278b32a342fSShri Abhyankar    option in the database.
127953acd3b1SBarry Smith 
12803f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith    Input Parameters:
128353acd3b1SBarry Smith +  opt - the option one is seeking
128453acd3b1SBarry Smith .  text - short string describing option
128553acd3b1SBarry Smith .  man - manual page for option
1286f8a50e2bSBarry Smith -  n - maximum number of values
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith    Output Parameter:
128953acd3b1SBarry Smith +  value - location to copy values
1290f8a50e2bSBarry Smith .  n - actual number of values found
129153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
129253acd3b1SBarry Smith 
129353acd3b1SBarry Smith    Level: beginner
129453acd3b1SBarry Smith 
129553acd3b1SBarry Smith    Notes:
1296b32a342fSShri Abhyankar    The array can be passed as
1297b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12980fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12990fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13000fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1301b32a342fSShri Abhyankar 
1302b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
130353acd3b1SBarry Smith 
130453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
130553acd3b1SBarry Smith 
1306b32a342fSShri Abhyankar    Concepts: options database^array of ints
130753acd3b1SBarry Smith 
130853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1309acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
131053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
131153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1312acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1313a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
131453acd3b1SBarry Smith @*/
13157087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
131653acd3b1SBarry Smith {
131753acd3b1SBarry Smith   PetscErrorCode ierr;
131853acd3b1SBarry Smith   PetscInt       i;
1319e26ddf31SBarry Smith   PetscOptions   amsopt;
132053acd3b1SBarry Smith 
132153acd3b1SBarry Smith   PetscFunctionBegin;
1322b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1323e26ddf31SBarry Smith     PetscInt *vals;
1324e26ddf31SBarry Smith 
1325e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1326e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1327e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1328e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1329e26ddf31SBarry Smith     amsopt->arraylength = *n;
1330e26ddf31SBarry Smith   }
133153acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
133261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
133353acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
133453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
133553acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
133653acd3b1SBarry Smith     }
13372aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133853acd3b1SBarry Smith   }
133953acd3b1SBarry Smith   PetscFunctionReturn(0);
134053acd3b1SBarry Smith }
134153acd3b1SBarry Smith 
134253acd3b1SBarry Smith #undef __FUNCT__
134353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
134453acd3b1SBarry Smith /*@C
134553acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
134653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
134753acd3b1SBarry Smith    no intervening spaces.
134853acd3b1SBarry Smith 
13493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
135053acd3b1SBarry Smith 
135153acd3b1SBarry Smith    Input Parameters:
135253acd3b1SBarry Smith +  opt - the option one is seeking
135353acd3b1SBarry Smith .  text - short string describing option
135453acd3b1SBarry Smith .  man - manual page for option
135553acd3b1SBarry Smith -  nmax - maximum number of strings
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith    Output Parameter:
135853acd3b1SBarry Smith +  value - location to copy strings
135953acd3b1SBarry Smith .  nmax - actual number of strings found
136053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
136153acd3b1SBarry Smith 
136253acd3b1SBarry Smith    Level: beginner
136353acd3b1SBarry Smith 
136453acd3b1SBarry Smith    Notes:
136553acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
136653acd3b1SBarry Smith    strings returned by this function.
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136953acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
137053acd3b1SBarry Smith 
137153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
137253acd3b1SBarry Smith 
137353acd3b1SBarry Smith    Concepts: options database^array of strings
137453acd3b1SBarry Smith 
137553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1376acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1379acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1380a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
138153acd3b1SBarry Smith @*/
13827087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
138353acd3b1SBarry Smith {
138453acd3b1SBarry Smith   PetscErrorCode ierr;
13851ae3d29cSBarry Smith   PetscOptions   amsopt;
138653acd3b1SBarry Smith 
138753acd3b1SBarry Smith   PetscFunctionBegin;
13881ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13891ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13901ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1391a297a907SKarl Rupp 
13921ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13931ae3d29cSBarry Smith   }
139453acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
139561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13962aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139753acd3b1SBarry Smith   }
139853acd3b1SBarry Smith   PetscFunctionReturn(0);
139953acd3b1SBarry Smith }
140053acd3b1SBarry Smith 
1401e2446a98SMatthew Knepley #undef __FUNCT__
1402acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1403e2446a98SMatthew Knepley /*@C
1404acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1405e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1406e2446a98SMatthew Knepley    no intervening spaces.
1407e2446a98SMatthew Knepley 
14083f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1409e2446a98SMatthew Knepley 
1410e2446a98SMatthew Knepley    Input Parameters:
1411e2446a98SMatthew Knepley +  opt - the option one is seeking
1412e2446a98SMatthew Knepley .  text - short string describing option
1413e2446a98SMatthew Knepley .  man - manual page for option
1414e2446a98SMatthew Knepley -  nmax - maximum number of values
1415e2446a98SMatthew Knepley 
1416e2446a98SMatthew Knepley    Output Parameter:
1417e2446a98SMatthew Knepley +  value - location to copy values
1418e2446a98SMatthew Knepley .  nmax - actual number of values found
1419e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1420e2446a98SMatthew Knepley 
1421e2446a98SMatthew Knepley    Level: beginner
1422e2446a98SMatthew Knepley 
1423e2446a98SMatthew Knepley    Notes:
1424e2446a98SMatthew Knepley    The user should pass in an array of doubles
1425e2446a98SMatthew Knepley 
1426e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1427e2446a98SMatthew Knepley 
1428e2446a98SMatthew Knepley    Concepts: options database^array of strings
1429e2446a98SMatthew Knepley 
1430e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1431acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1432e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1433e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1434acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1435a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1436e2446a98SMatthew Knepley @*/
14377087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1438e2446a98SMatthew Knepley {
1439e2446a98SMatthew Knepley   PetscErrorCode ierr;
1440e2446a98SMatthew Knepley   PetscInt       i;
14411ae3d29cSBarry Smith   PetscOptions   amsopt;
1442e2446a98SMatthew Knepley 
1443e2446a98SMatthew Knepley   PetscFunctionBegin;
14441ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1445ace3abfcSBarry Smith     PetscBool *vals;
14461ae3d29cSBarry Smith 
14477781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1448ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1449ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14501ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14511ae3d29cSBarry Smith     amsopt->arraylength = *n;
14521ae3d29cSBarry Smith   }
1453acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1454e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1455e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1456e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1457e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1458e2446a98SMatthew Knepley     }
14592aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1460e2446a98SMatthew Knepley   }
1461e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1462e2446a98SMatthew Knepley }
1463e2446a98SMatthew Knepley 
14648cc676e6SMatthew G Knepley #undef __FUNCT__
14658cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14668cc676e6SMatthew G Knepley /*@C
14678cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14688cc676e6SMatthew G Knepley 
14698cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14708cc676e6SMatthew G Knepley 
14718cc676e6SMatthew G Knepley    Input Parameters:
14728cc676e6SMatthew G Knepley +  opt - option name
14738cc676e6SMatthew G Knepley .  text - short string that describes the option
14748cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14758cc676e6SMatthew G Knepley 
14768cc676e6SMatthew G Knepley    Output Parameter:
14778cc676e6SMatthew G Knepley +  viewer - the viewer
14788cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14798cc676e6SMatthew G Knepley 
14808cc676e6SMatthew G Knepley    Level: beginner
14818cc676e6SMatthew G Knepley 
14828cc676e6SMatthew G Knepley    Concepts: options database^has int
14838cc676e6SMatthew G Knepley 
14848cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14858cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14868cc676e6SMatthew 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
14878cc676e6SMatthew G Knepley $                                     about the object to standard out
14888cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14898cc676e6SMatthew G Knepley $       draw
14908cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14918cc676e6SMatthew G Knepley 
1492cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14938cc676e6SMatthew G Knepley 
14948cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14958cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14968cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14978cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14988cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14998cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1500a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15018cc676e6SMatthew G Knepley @*/
1502cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15038cc676e6SMatthew G Knepley {
15048cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15058cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15068cc676e6SMatthew G Knepley 
15078cc676e6SMatthew G Knepley   PetscFunctionBegin;
15088cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15098cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
1510*64facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
1511c979a496SBarry Smith     amsopt->data = (void*)strdup("");
15128cc676e6SMatthew G Knepley   }
1513cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15148cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15158cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15168cc676e6SMatthew G Knepley   }
15178cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15188cc676e6SMatthew G Knepley }
15198cc676e6SMatthew G Knepley 
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith #undef __FUNCT__
152253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
152353acd3b1SBarry Smith /*@C
1524b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
152553acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
152653acd3b1SBarry Smith 
15273f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
152853acd3b1SBarry Smith 
152953acd3b1SBarry Smith    Input Parameter:
153053acd3b1SBarry Smith .   head - the heading text
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith 
153353acd3b1SBarry Smith    Level: intermediate
153453acd3b1SBarry Smith 
153553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
153653acd3b1SBarry Smith 
1537b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
153853acd3b1SBarry Smith 
153953acd3b1SBarry Smith    Concepts: options database^subheading
154053acd3b1SBarry Smith 
154153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1542acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
154353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
154453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1545acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1546a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
154753acd3b1SBarry Smith @*/
15487087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154953acd3b1SBarry Smith {
155053acd3b1SBarry Smith   PetscErrorCode ierr;
155153acd3b1SBarry Smith 
155253acd3b1SBarry Smith   PetscFunctionBegin;
155361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
155453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
155553acd3b1SBarry Smith   }
155653acd3b1SBarry Smith   PetscFunctionReturn(0);
155753acd3b1SBarry Smith }
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith 
156053acd3b1SBarry Smith 
156153acd3b1SBarry Smith 
156253acd3b1SBarry Smith 
156353acd3b1SBarry Smith 
1564