xref: /petsc/src/sys/objects/aoptions.c (revision afcb2eb581850588619400ff311bed89835ec555)
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 
8*afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
953acd3b1SBarry Smith 
102aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
112aa6d131SJed Brown 
1253acd3b1SBarry Smith /*
1353acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
143fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1553acd3b1SBarry Smith 
1653acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1753acd3b1SBarry Smith */
18f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
1953acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2053acd3b1SBarry Smith 
2153acd3b1SBarry Smith #undef __FUNCT__
2253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2353acd3b1SBarry Smith /*
2453acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2553acd3b1SBarry Smith */
2653acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2753acd3b1SBarry Smith {
2853acd3b1SBarry Smith   PetscErrorCode ierr;
2953acd3b1SBarry Smith 
3053acd3b1SBarry Smith   PetscFunctionBegin;
3153acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3253acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
336356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
34a297a907SKarl Rupp 
35c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3653acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
37c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3853acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith 
400298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4153acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4261b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4353acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4453acd3b1SBarry Smith     }
4561b37b28SSatish Balay   }
4653acd3b1SBarry Smith   PetscFunctionReturn(0);
4753acd3b1SBarry Smith }
4853acd3b1SBarry Smith 
493194b578SJed Brown #undef __FUNCT__
503194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
513194b578SJed Brown /*
523194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
533194b578SJed Brown */
543194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
553194b578SJed Brown {
563194b578SJed Brown   PetscErrorCode ierr;
573194b578SJed Brown   char           title[256];
583194b578SJed Brown   PetscBool      flg;
593194b578SJed Brown 
603194b578SJed Brown   PetscFunctionBegin;
613194b578SJed Brown   PetscValidHeader(obj,1);
623194b578SJed Brown   PetscOptionsObject.object         = obj;
633194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
64a297a907SKarl Rupp 
653194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
663194b578SJed Brown   if (flg) {
678caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
683194b578SJed Brown   } else {
698caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
703194b578SJed Brown   }
713194b578SJed Brown   ierr = PetscOptionsBegin_Private(obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
723194b578SJed Brown   PetscFunctionReturn(0);
733194b578SJed Brown }
743194b578SJed Brown 
7553acd3b1SBarry Smith /*
7653acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7753acd3b1SBarry Smith */
7853acd3b1SBarry Smith #undef __FUNCT__
7953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
80e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8153acd3b1SBarry Smith {
8253acd3b1SBarry Smith   int          ierr;
8353acd3b1SBarry Smith   PetscOptions next;
843be6e4c3SJed Brown   PetscBool    valid;
8553acd3b1SBarry Smith 
8653acd3b1SBarry Smith   PetscFunctionBegin;
873be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
883be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
893be6e4c3SJed Brown 
906e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9153acd3b1SBarry Smith   (*amsopt)->next = 0;
9253acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
936356e834SBarry Smith   (*amsopt)->type = t;
9453acd3b1SBarry Smith   (*amsopt)->data = 0;
9561b37b28SSatish Balay 
9653acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
986356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
9953acd3b1SBarry Smith 
100a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
101a297a907SKarl Rupp   else {
10253acd3b1SBarry Smith     next = PetscOptionsObject.next;
10353acd3b1SBarry Smith     while (next->next) next = next->next;
10453acd3b1SBarry Smith     next->next = *amsopt;
10553acd3b1SBarry Smith   }
10653acd3b1SBarry Smith   PetscFunctionReturn(0);
10753acd3b1SBarry Smith }
10853acd3b1SBarry Smith 
10953acd3b1SBarry Smith #undef __FUNCT__
110aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
111aee2cecaSBarry Smith /*
1123fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1133fc1eb6aSBarry Smith 
1143fc1eb6aSBarry Smith     Collective on MPI_Comm
1153fc1eb6aSBarry Smith 
1163fc1eb6aSBarry Smith    Input Parameters:
1173fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1183fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1193fc1eb6aSBarry Smith -     str - location to store input
120aee2cecaSBarry Smith 
121aee2cecaSBarry Smith     Bugs:
122aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
123aee2cecaSBarry Smith 
124aee2cecaSBarry Smith */
1253fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
126aee2cecaSBarry Smith {
127330cf3c9SBarry Smith   size_t         i;
128aee2cecaSBarry Smith   char           c;
1293fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
130aee2cecaSBarry Smith   PetscErrorCode ierr;
131aee2cecaSBarry Smith 
132aee2cecaSBarry Smith   PetscFunctionBegin;
133aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
134aee2cecaSBarry Smith   if (!rank) {
135aee2cecaSBarry Smith     c = (char) getchar();
136aee2cecaSBarry Smith     i = 0;
137aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
138aee2cecaSBarry Smith       str[i++] = c;
139aee2cecaSBarry Smith       c = (char)getchar();
140aee2cecaSBarry Smith     }
141aee2cecaSBarry Smith     str[i] = 0;
142aee2cecaSBarry Smith   }
1434dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1443fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
145aee2cecaSBarry Smith   PetscFunctionReturn(0);
146aee2cecaSBarry Smith }
147aee2cecaSBarry Smith 
148aee2cecaSBarry Smith #undef __FUNCT__
149aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
150aee2cecaSBarry Smith /*
1513cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
152aee2cecaSBarry Smith 
153aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
154aee2cecaSBarry Smith 
155aee2cecaSBarry Smith     Bugs:
156aee2cecaSBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
1573cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
158aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
159aee2cecaSBarry Smith 
1603cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1613cc1e11dSBarry Smith      address space and communicating with the PETSc program
1623cc1e11dSBarry Smith 
163aee2cecaSBarry Smith */
164aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1656356e834SBarry Smith {
1666356e834SBarry Smith   PetscErrorCode ierr;
1676356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1686356e834SBarry Smith   char           str[512];
169a4404d99SBarry Smith   PetscInt       id;
170a4404d99SBarry Smith   PetscReal      ir,*valr;
171330cf3c9SBarry Smith   PetscInt       *vald;
172330cf3c9SBarry Smith   size_t         i;
1736356e834SBarry Smith 
174e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1756356e834SBarry Smith   while (next) {
1766356e834SBarry Smith     switch (next->type) {
1776356e834SBarry Smith     case OPTION_HEAD:
1786356e834SBarry Smith       break;
179e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
180e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
181e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
182e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
183e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
184e26ddf31SBarry Smith         if (i < next->arraylength-1) {
185e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
186e26ddf31SBarry Smith         }
187e26ddf31SBarry Smith       }
188e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
189e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
190e26ddf31SBarry Smith       if (str[0]) {
191e26ddf31SBarry Smith         PetscToken token;
192e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
193e26ddf31SBarry Smith         size_t     len;
194e26ddf31SBarry Smith         char       *value;
195ace3abfcSBarry Smith         PetscBool  foundrange;
196e26ddf31SBarry Smith 
197e26ddf31SBarry Smith         next->set = PETSC_TRUE;
198e26ddf31SBarry Smith         value     = str;
199e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
200e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
201e26ddf31SBarry Smith         while (n < nmax) {
202e26ddf31SBarry Smith           if (!value) break;
203e26ddf31SBarry Smith 
204e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
205e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
206e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
207e26ddf31SBarry Smith           if (value[0] == '-') i=2;
208e26ddf31SBarry Smith           else i=1;
209330cf3c9SBarry Smith           for (;i<len; i++) {
210e26ddf31SBarry Smith             if (value[i] == '-') {
211e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
212e26ddf31SBarry Smith               value[i] = 0;
213cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
214cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
215e32f2f54SBarry 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);
216e32f2f54SBarry 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);
217e26ddf31SBarry Smith               for (; start<end; start++) {
218e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
219e26ddf31SBarry Smith               }
220e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
221e26ddf31SBarry Smith               break;
222e26ddf31SBarry Smith             }
223e26ddf31SBarry Smith           }
224e26ddf31SBarry Smith           if (!foundrange) {
225cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
226e26ddf31SBarry Smith             dvalue++;
227e26ddf31SBarry Smith             n++;
228e26ddf31SBarry Smith           }
229e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
230e26ddf31SBarry Smith         }
2318c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
232e26ddf31SBarry Smith       }
233e26ddf31SBarry Smith       break;
234e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
235e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
236e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
237e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
238e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
239e26ddf31SBarry Smith         if (i < next->arraylength-1) {
240e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
241e26ddf31SBarry Smith         }
242e26ddf31SBarry Smith       }
243e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
244e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
245e26ddf31SBarry Smith       if (str[0]) {
246e26ddf31SBarry Smith         PetscToken token;
247e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
248e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
249e26ddf31SBarry Smith         char       *value;
250e26ddf31SBarry Smith 
251e26ddf31SBarry Smith         next->set = PETSC_TRUE;
252e26ddf31SBarry Smith         value     = str;
253e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
254e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
255e26ddf31SBarry Smith         while (n < nmax) {
256e26ddf31SBarry Smith           if (!value) break;
257cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
258e26ddf31SBarry Smith           dvalue++;
259e26ddf31SBarry Smith           n++;
260e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
261e26ddf31SBarry Smith         }
2628c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
263e26ddf31SBarry Smith       }
264e26ddf31SBarry Smith       break;
2656356e834SBarry Smith     case OPTION_INT:
266e26ddf31SBarry 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);
2673fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2683fc1eb6aSBarry Smith       if (str[0]) {
269c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
270c272547aSJed Brown         sscanf(str,"%lld",&id);
271c272547aSJed Brown #else
272aee2cecaSBarry Smith         sscanf(str,"%d",&id);
273c272547aSJed Brown #endif
274aee2cecaSBarry Smith         next->set = PETSC_TRUE;
275a297a907SKarl Rupp 
276aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
277aee2cecaSBarry Smith       }
278aee2cecaSBarry Smith       break;
279aee2cecaSBarry Smith     case OPTION_REAL:
280e26ddf31SBarry 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);
2813fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2823fc1eb6aSBarry Smith       if (str[0]) {
283ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
284a4404d99SBarry Smith         sscanf(str,"%e",&ir);
285ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
286aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
287ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
288d9822059SBarry Smith         ir = strtoflt128(str,0);
289d9822059SBarry Smith #else
290513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
291a4404d99SBarry Smith #endif
292aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
293aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
294aee2cecaSBarry Smith       }
295aee2cecaSBarry Smith       break;
296aee2cecaSBarry Smith     case OPTION_LOGICAL:
297aee2cecaSBarry Smith     case OPTION_STRING:
298e26ddf31SBarry 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);
2993fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3003fc1eb6aSBarry Smith       if (str[0]) {
301aee2cecaSBarry Smith         next->set = PETSC_TRUE;
302a297a907SKarl Rupp 
303330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3046356e834SBarry Smith       }
3056356e834SBarry Smith       break;
3063cc1e11dSBarry Smith     case OPTION_LIST:
307140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3083cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3093cc1e11dSBarry Smith       if (str[0]) {
3103cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3113cc1e11dSBarry Smith         next->set = PETSC_TRUE;
312330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3133cc1e11dSBarry Smith       }
3143cc1e11dSBarry Smith       break;
315b432afdaSMatthew Knepley     default:
316b432afdaSMatthew Knepley       break;
3176356e834SBarry Smith     }
3186356e834SBarry Smith     next = next->next;
3196356e834SBarry Smith   }
3206356e834SBarry Smith   PetscFunctionReturn(0);
3216356e834SBarry Smith }
3226356e834SBarry Smith 
323b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
3241ae3d29cSBarry Smith #define CHKERRAMS(err)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"AMS Error: %s",msg);}
3251ae3d29cSBarry Smith #define CHKERRAMSFieldName(err,fn)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Fieldname %s, AMS Error: %s",fn,msg);}
326d5649816SBarry Smith 
327d5649816SBarry Smith static int count = 0;
328d5649816SBarry Smith 
329b3506946SBarry Smith #undef __FUNCT__
330d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSDestroy"
331d5649816SBarry Smith PetscErrorCode PetscOptionsAMSDestroy(void)
332d5649816SBarry Smith {
333d5649816SBarry Smith   PetscErrorCode ierr;
334d5649816SBarry Smith   AMS_Comm       acomm = -1;
335d5649816SBarry Smith   AMS_Memory     amem  = -1;
336d5649816SBarry Smith   char           options[16];
337d5649816SBarry Smith   const char     *string = "Exit";
338d5649816SBarry Smith 
339d5649816SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
340d5649816SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
341d5649816SBarry Smith   sprintf(options,"Options_%d",count++);
342d5649816SBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
343d5649816SBarry Smith   ierr = AMS_Memory_add_field(amem,"Exit",&string,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"Exit");
344d5649816SBarry Smith 
345d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
346d5649816SBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
347d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
348d5649816SBarry Smith   /* wait until accessor has unlocked the memory */
349d5649816SBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
350d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
351d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
352d5649816SBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
353d5649816SBarry Smith   PetscFunctionReturn(0);
354d5649816SBarry Smith }
355d5649816SBarry Smith 
356d5649816SBarry Smith #undef __FUNCT__
357d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSInput"
358b3506946SBarry Smith /*
359d5649816SBarry Smith     PetscOptionsAMSInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the AMS
360b3506946SBarry Smith 
361b3506946SBarry Smith     Bugs:
362b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
363b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
364b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
365b3506946SBarry Smith 
366b3506946SBarry Smith 
367b3506946SBarry Smith */
368d5649816SBarry Smith PetscErrorCode PetscOptionsAMSInput()
369b3506946SBarry Smith {
370b3506946SBarry Smith   PetscErrorCode ierr;
371b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
372d5649816SBarry Smith   static int     mancount = 0;
373b3506946SBarry Smith   char           options[16];
374b3506946SBarry Smith   AMS_Comm       acomm         = -1;
375b3506946SBarry Smith   AMS_Memory     amem          = -1;
376ace3abfcSBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
3779f32e415SBarry Smith   char           manname[16];
378b3506946SBarry Smith 
379b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
380b3506946SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
381b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
3821ae3d29cSBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
3831ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
384a297a907SKarl Rupp 
3851bc75a8dSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* AMS will change this, so cannot pass prefix directly */
3861bc75a8dSBarry Smith 
3871ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,PetscOptionsObject.title,&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,PetscOptionsObject.title);
3880c74a584SJed Brown   /* ierr = AMS_Memory_add_field(amem,"mansec",&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMS(ierr); */
3891ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,"ChangedMethod",&changedmethod,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"ChangedMethod");
390b3506946SBarry Smith 
391b3506946SBarry Smith   while (next) {
3921ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,next->option,&next->set,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->option);
3931bc75a8dSBarry Smith     ierr =  PetscMalloc(sizeof(char*),&next->pman);CHKERRQ(ierr);
394a297a907SKarl Rupp 
3951bc75a8dSBarry Smith     *(char**)next->pman = next->man;
3969f32e415SBarry Smith     sprintf(manname,"man_%d",mancount++);
3971ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,manname,next->pman,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,manname);
3989f32e415SBarry Smith 
399b3506946SBarry Smith     switch (next->type) {
400b3506946SBarry Smith     case OPTION_HEAD:
401b3506946SBarry Smith       break;
402b3506946SBarry Smith     case OPTION_INT_ARRAY:
4031ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
404b3506946SBarry Smith       break;
405b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4061ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
407b3506946SBarry Smith       break;
408b3506946SBarry Smith     case OPTION_INT:
4091ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
410b3506946SBarry Smith       break;
411b3506946SBarry Smith     case OPTION_REAL:
4121ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
413b3506946SBarry Smith       break;
414b3506946SBarry Smith     case OPTION_LOGICAL:
4151ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4161ae3d29cSBarry Smith       break;
4171ae3d29cSBarry Smith     case OPTION_LOGICAL_ARRAY:
4181ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
41971f08665SBarry Smith       break;
420b3506946SBarry Smith     case OPTION_STRING:
4211ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4221ae3d29cSBarry Smith       break;
4231ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4241ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
425b3506946SBarry Smith       break;
426b3506946SBarry Smith     case OPTION_LIST:
42771f08665SBarry Smith       {PetscInt ntext;
42871f08665SBarry Smith       char      ldefault[128];
42971f08665SBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
43071f08665SBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4311ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
432140e18c1SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
433d5649816SBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext-1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4341ae3d29cSBarry Smith       break;}
4351ae3d29cSBarry Smith     case OPTION_ELIST:
436d5649816SBarry Smith       {PetscInt ntext = next->nlist;
4371ae3d29cSBarry Smith       char      ldefault[128];
4381ae3d29cSBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
4391ae3d29cSBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4401ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
4411ae3d29cSBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4421ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4431ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
44471f08665SBarry Smith       break;}
445b3506946SBarry Smith     default:
446b3506946SBarry Smith       break;
447b3506946SBarry Smith     }
448b3506946SBarry Smith     next = next->next;
449b3506946SBarry Smith   }
450b3506946SBarry Smith 
4511ae3d29cSBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
4521ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
453b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4541ae3d29cSBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
4551ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
456b3506946SBarry Smith 
457b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
458b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
459b3506946SBarry Smith 
4601ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
4611ae3d29cSBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
462b3506946SBarry Smith   PetscFunctionReturn(0);
463b3506946SBarry Smith }
464b3506946SBarry Smith #endif
465b3506946SBarry Smith 
4666356e834SBarry Smith #undef __FUNCT__
46753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
46853acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
46953acd3b1SBarry Smith {
47053acd3b1SBarry Smith   PetscErrorCode ierr;
4716356e834SBarry Smith   PetscOptions   last;
4726356e834SBarry Smith   char           option[256],value[1024],tmp[32];
473330cf3c9SBarry Smith   size_t         j;
47453acd3b1SBarry Smith 
47553acd3b1SBarry Smith   PetscFunctionBegin;
4761bc75a8dSBarry Smith   CHKMEMQ;
477aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
478b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
479b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
480d5649816SBarry Smith       ierr = PetscOptionsAMSInput();CHKERRQ(ierr);
481b3506946SBarry Smith #else
48271f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
483b3506946SBarry Smith #endif
484aee2cecaSBarry Smith     }
485aee2cecaSBarry Smith   }
4866356e834SBarry Smith 
487c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
488c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4896356e834SBarry Smith 
4906356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4916356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49261b37b28SSatish Balay   /* reset alreadyprinted flag */
49361b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4943194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4950298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4966356e834SBarry Smith 
4976356e834SBarry Smith   while (PetscOptionsObject.next) {
4986356e834SBarry Smith     if (PetscOptionsObject.next->set) {
4996356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5006356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5016356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5026356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5036356e834SBarry Smith       } else {
5046356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5056356e834SBarry Smith       }
5066356e834SBarry Smith 
5076356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5086356e834SBarry Smith       case OPTION_HEAD:
5096356e834SBarry Smith         break;
510e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
511e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
512e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
513e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
514e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
515e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
516e26ddf31SBarry Smith         }
517e26ddf31SBarry Smith         break;
5186356e834SBarry Smith       case OPTION_INT:
5197a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5206356e834SBarry Smith         break;
5216356e834SBarry Smith       case OPTION_REAL:
5227a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5236356e834SBarry Smith         break;
5246356e834SBarry Smith       case OPTION_REAL_ARRAY:
5257a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5266356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5277a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5286356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5296356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5306356e834SBarry Smith         }
5316356e834SBarry Smith         break;
5326356e834SBarry Smith       case OPTION_LOGICAL:
53371f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5346356e834SBarry Smith         break;
5351ae3d29cSBarry Smith       case OPTION_LOGICAL_ARRAY:
536ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5371ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
538ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5391ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5401ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5411ae3d29cSBarry Smith         }
5421ae3d29cSBarry Smith         break;
5436356e834SBarry Smith       case OPTION_LIST:
5441ae3d29cSBarry Smith       case OPTION_ELIST:
54571f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5466356e834SBarry Smith         break;
5471ae3d29cSBarry Smith       case OPTION_STRING:
54871f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5491ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5501ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5511ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5521ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5531ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5541ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5551ae3d29cSBarry Smith         }
5566356e834SBarry Smith         break;
5576356e834SBarry Smith       }
5586356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5596356e834SBarry Smith     }
560503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
561503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5626356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56305b42c5fSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
56471f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
565a297a907SKarl Rupp 
5666356e834SBarry Smith     last                    = PetscOptionsObject.next;
5676356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5686356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5691bc75a8dSBarry Smith     CHKMEMQ;
5706356e834SBarry Smith   }
5711bc75a8dSBarry Smith   CHKMEMQ;
5726356e834SBarry Smith   PetscOptionsObject.next = 0;
57353acd3b1SBarry Smith   PetscFunctionReturn(0);
57453acd3b1SBarry Smith }
57553acd3b1SBarry Smith 
57653acd3b1SBarry Smith #undef __FUNCT__
57753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
57853acd3b1SBarry Smith /*@C
57953acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58053acd3b1SBarry Smith 
5813f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58253acd3b1SBarry Smith 
58353acd3b1SBarry Smith    Input Parameters:
58453acd3b1SBarry Smith +  opt - option name
58553acd3b1SBarry Smith .  text - short string that describes the option
58653acd3b1SBarry Smith .  man - manual page with additional information on option
58753acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
58853acd3b1SBarry Smith -  defaultv - the default (current) value
58953acd3b1SBarry Smith 
59053acd3b1SBarry Smith    Output Parameter:
59153acd3b1SBarry Smith +  value - the  value to return
59253acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
59353acd3b1SBarry Smith 
59453acd3b1SBarry Smith    Level: beginner
59553acd3b1SBarry Smith 
59653acd3b1SBarry Smith    Concepts: options database
59753acd3b1SBarry Smith 
59853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
59953acd3b1SBarry Smith 
60053acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60153acd3b1SBarry Smith 
60253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
603acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
604acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
607acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
60853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
60953acd3b1SBarry Smith @*/
6107087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61153acd3b1SBarry Smith {
61253acd3b1SBarry Smith   PetscErrorCode ierr;
61353acd3b1SBarry Smith   PetscInt       ntext = 0;
614aa5bb8c0SSatish Balay   PetscInt       tval;
615ace3abfcSBarry Smith   PetscBool      tflg;
61653acd3b1SBarry Smith 
61753acd3b1SBarry Smith   PetscFunctionBegin;
61853acd3b1SBarry Smith   while (list[ntext++]) {
619e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62053acd3b1SBarry Smith   }
621e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62253acd3b1SBarry Smith   ntext -= 3;
623aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
624aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
625aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
626aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
62753acd3b1SBarry Smith   PetscFunctionReturn(0);
62853acd3b1SBarry Smith }
62953acd3b1SBarry Smith 
63053acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63153acd3b1SBarry Smith #undef __FUNCT__
63253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63353acd3b1SBarry Smith /*@C
63453acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63553acd3b1SBarry Smith 
6363f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63753acd3b1SBarry Smith 
63853acd3b1SBarry Smith    Input Parameters:
63953acd3b1SBarry Smith +  opt - option name
64053acd3b1SBarry Smith .  text - short string that describes the option
64153acd3b1SBarry Smith .  man - manual page with additional information on option
64253acd3b1SBarry Smith -  defaultv - the default (current) value
64353acd3b1SBarry Smith 
64453acd3b1SBarry Smith    Output Parameter:
64553acd3b1SBarry Smith +  value - the integer value to return
64653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
64753acd3b1SBarry Smith 
64853acd3b1SBarry Smith    Level: beginner
64953acd3b1SBarry Smith 
65053acd3b1SBarry Smith    Concepts: options database^has int
65153acd3b1SBarry Smith 
65253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65353acd3b1SBarry Smith 
65453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
655acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
656acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
659acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
66053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
66153acd3b1SBarry Smith @*/
6627087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66353acd3b1SBarry Smith {
66453acd3b1SBarry Smith   PetscErrorCode ierr;
6656356e834SBarry Smith   PetscOptions   amsopt;
66653acd3b1SBarry Smith 
66753acd3b1SBarry Smith   PetscFunctionBegin;
668b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6696356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6706356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
671a297a907SKarl Rupp 
6726356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
673af6d86caSBarry Smith   }
67453acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6762aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
67753acd3b1SBarry Smith   }
67853acd3b1SBarry Smith   PetscFunctionReturn(0);
67953acd3b1SBarry Smith }
68053acd3b1SBarry Smith 
68153acd3b1SBarry Smith #undef __FUNCT__
68253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68353acd3b1SBarry Smith /*@C
68453acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68553acd3b1SBarry Smith 
6863f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith    Input Parameters:
68953acd3b1SBarry Smith +  opt - option name
69053acd3b1SBarry Smith .  text - short string that describes the option
69153acd3b1SBarry Smith .  man - manual page with additional information on option
692bcbf2dc5SJed Brown .  defaultv - the default (current) value
693bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69453acd3b1SBarry Smith 
69553acd3b1SBarry Smith    Output Parameter:
69653acd3b1SBarry Smith +  value - the value to return
69753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69853acd3b1SBarry Smith 
69953acd3b1SBarry Smith    Level: beginner
70053acd3b1SBarry Smith 
70153acd3b1SBarry Smith    Concepts: options database^has int
70253acd3b1SBarry Smith 
70353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70453acd3b1SBarry Smith 
7057fccdfe4SBarry Smith    Even if the user provided no string (for example -optionname -someotheroption) the flag is set to PETSC_TRUE (and the string is fulled with nulls).
7067fccdfe4SBarry Smith 
70753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
708acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
709acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
712acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71453acd3b1SBarry Smith @*/
7157087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71653acd3b1SBarry Smith {
71753acd3b1SBarry Smith   PetscErrorCode ierr;
718aee2cecaSBarry Smith   PetscOptions   amsopt;
71953acd3b1SBarry Smith 
72053acd3b1SBarry Smith   PetscFunctionBegin;
721b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
722aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
72371f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
724a297a907SKarl Rupp 
72571f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
726af6d86caSBarry Smith   }
72753acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
72861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7292aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73053acd3b1SBarry Smith   }
73153acd3b1SBarry Smith   PetscFunctionReturn(0);
73253acd3b1SBarry Smith }
73353acd3b1SBarry Smith 
73453acd3b1SBarry Smith #undef __FUNCT__
73553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73653acd3b1SBarry Smith /*@C
73753acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
73853acd3b1SBarry Smith 
7393f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74053acd3b1SBarry Smith 
74153acd3b1SBarry Smith    Input Parameters:
74253acd3b1SBarry Smith +  opt - option name
74353acd3b1SBarry Smith .  text - short string that describes the option
74453acd3b1SBarry Smith .  man - manual page with additional information on option
74553acd3b1SBarry Smith -  defaultv - the default (current) value
74653acd3b1SBarry Smith 
74753acd3b1SBarry Smith    Output Parameter:
74853acd3b1SBarry Smith +  value - the value to return
74953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75053acd3b1SBarry Smith 
75153acd3b1SBarry Smith    Level: beginner
75253acd3b1SBarry Smith 
75353acd3b1SBarry Smith    Concepts: options database^has int
75453acd3b1SBarry Smith 
75553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75653acd3b1SBarry Smith 
75753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
758acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
759acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
76053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
762acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
76353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
76453acd3b1SBarry Smith @*/
7657087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76653acd3b1SBarry Smith {
76753acd3b1SBarry Smith   PetscErrorCode ierr;
768538aa990SBarry Smith   PetscOptions   amsopt;
76953acd3b1SBarry Smith 
77053acd3b1SBarry Smith   PetscFunctionBegin;
771b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
772538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
773538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
774a297a907SKarl Rupp 
775538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
776538aa990SBarry Smith   }
77753acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
77861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7792aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78053acd3b1SBarry Smith   }
78153acd3b1SBarry Smith   PetscFunctionReturn(0);
78253acd3b1SBarry Smith }
78353acd3b1SBarry Smith 
78453acd3b1SBarry Smith #undef __FUNCT__
78553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78653acd3b1SBarry Smith /*@C
78753acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
78853acd3b1SBarry Smith 
7893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79053acd3b1SBarry Smith 
79153acd3b1SBarry Smith    Input Parameters:
79253acd3b1SBarry Smith +  opt - option name
79353acd3b1SBarry Smith .  text - short string that describes the option
79453acd3b1SBarry Smith .  man - manual page with additional information on option
79553acd3b1SBarry Smith -  defaultv - the default (current) value
79653acd3b1SBarry Smith 
79753acd3b1SBarry Smith    Output Parameter:
79853acd3b1SBarry Smith +  value - the value to return
79953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80053acd3b1SBarry Smith 
80153acd3b1SBarry Smith    Level: beginner
80253acd3b1SBarry Smith 
80353acd3b1SBarry Smith    Concepts: options database^has int
80453acd3b1SBarry Smith 
80553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80653acd3b1SBarry Smith 
80753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
808acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
809acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
812acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
81353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
81453acd3b1SBarry Smith @*/
8157087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81653acd3b1SBarry Smith {
81753acd3b1SBarry Smith   PetscErrorCode ierr;
81853acd3b1SBarry Smith 
81953acd3b1SBarry Smith   PetscFunctionBegin;
82053acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
82153acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82253acd3b1SBarry Smith #else
82353acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82453acd3b1SBarry Smith #endif
82553acd3b1SBarry Smith   PetscFunctionReturn(0);
82653acd3b1SBarry Smith }
82753acd3b1SBarry Smith 
82853acd3b1SBarry Smith #undef __FUNCT__
82953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
83053acd3b1SBarry Smith /*@C
83190d69ab7SBarry 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
83290d69ab7SBarry Smith                       its value is set to false.
83353acd3b1SBarry Smith 
8343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith    Input Parameters:
83753acd3b1SBarry Smith +  opt - option name
83853acd3b1SBarry Smith .  text - short string that describes the option
83953acd3b1SBarry Smith -  man - manual page with additional information on option
84053acd3b1SBarry Smith 
84153acd3b1SBarry Smith    Output Parameter:
84253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
84353acd3b1SBarry Smith 
84453acd3b1SBarry Smith    Level: beginner
84553acd3b1SBarry Smith 
84653acd3b1SBarry Smith    Concepts: options database^has int
84753acd3b1SBarry Smith 
84853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
84953acd3b1SBarry Smith 
85053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
851acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
852acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
855acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
85653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
85753acd3b1SBarry Smith @*/
8587087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
85953acd3b1SBarry Smith {
86053acd3b1SBarry Smith   PetscErrorCode ierr;
8611ae3d29cSBarry Smith   PetscOptions   amsopt;
86253acd3b1SBarry Smith 
86353acd3b1SBarry Smith   PetscFunctionBegin;
8641ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8651ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
866ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
867a297a907SKarl Rupp 
868ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8691ae3d29cSBarry Smith   }
87053acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
87161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8722aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
87353acd3b1SBarry Smith   }
87453acd3b1SBarry Smith   PetscFunctionReturn(0);
87553acd3b1SBarry Smith }
87653acd3b1SBarry Smith 
87753acd3b1SBarry Smith #undef __FUNCT__
87853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
87953acd3b1SBarry Smith /*@C
88053acd3b1SBarry Smith      PetscOptionsList - Puts a list of option values that a single one may be selected from
88153acd3b1SBarry Smith 
8823f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88353acd3b1SBarry Smith 
88453acd3b1SBarry Smith    Input Parameters:
88553acd3b1SBarry Smith +  opt - option name
88653acd3b1SBarry Smith .  text - short string that describes the option
88753acd3b1SBarry Smith .  man - manual page with additional information on option
88853acd3b1SBarry Smith .  list - the possible choices
8893cc1e11dSBarry Smith .  defaultv - the default (current) value
8903cc1e11dSBarry Smith -  len - the length of the character array value
89153acd3b1SBarry Smith 
89253acd3b1SBarry Smith    Output Parameter:
89353acd3b1SBarry Smith +  value - the value to return
89453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith    Level: intermediate
89753acd3b1SBarry Smith 
89853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89953acd3b1SBarry Smith 
90053acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
90153acd3b1SBarry Smith 
90253acd3b1SBarry Smith    To get a listing of all currently specified options,
90388c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
90453acd3b1SBarry Smith 
90553acd3b1SBarry Smith    Concepts: options database^list
90653acd3b1SBarry Smith 
90753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
908acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
90953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
911acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
912171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
91353acd3b1SBarry Smith @*/
914140e18c1SBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91553acd3b1SBarry Smith {
91653acd3b1SBarry Smith   PetscErrorCode ierr;
9173cc1e11dSBarry Smith   PetscOptions   amsopt;
91853acd3b1SBarry Smith 
91953acd3b1SBarry Smith   PetscFunctionBegin;
920b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9213cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
92271f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
923a297a907SKarl Rupp 
92471f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
925a297a907SKarl Rupp 
9263cc1e11dSBarry Smith     amsopt->flist = list;
9273cc1e11dSBarry Smith   }
92853acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
92961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
930140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
93153acd3b1SBarry Smith   }
93253acd3b1SBarry Smith   PetscFunctionReturn(0);
93353acd3b1SBarry Smith }
93453acd3b1SBarry Smith 
93553acd3b1SBarry Smith #undef __FUNCT__
93653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
93753acd3b1SBarry Smith /*@C
93853acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
93953acd3b1SBarry Smith 
9403f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith    Input Parameters:
94353acd3b1SBarry Smith +  opt - option name
94453acd3b1SBarry Smith .  ltext - short string that describes the option
94553acd3b1SBarry Smith .  man - manual page with additional information on option
94653acd3b1SBarry Smith .  list - the possible choices
94753acd3b1SBarry Smith .  ntext - number of choices
94853acd3b1SBarry Smith -  defaultv - the default (current) value
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    Output Parameter:
95153acd3b1SBarry Smith +  value - the index of the value to return
95253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95353acd3b1SBarry Smith 
95453acd3b1SBarry Smith    Level: intermediate
95553acd3b1SBarry Smith 
95653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95753acd3b1SBarry Smith 
958140e18c1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFunctionList()
95953acd3b1SBarry Smith 
96053acd3b1SBarry Smith    Concepts: options database^list
96153acd3b1SBarry Smith 
96253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
963acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
966acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
967171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
96853acd3b1SBarry Smith @*/
9697087cfbeSBarry 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)
97053acd3b1SBarry Smith {
97153acd3b1SBarry Smith   PetscErrorCode ierr;
97253acd3b1SBarry Smith   PetscInt       i;
9731ae3d29cSBarry Smith   PetscOptions   amsopt;
97453acd3b1SBarry Smith 
97553acd3b1SBarry Smith   PetscFunctionBegin;
9761ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
977d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9781ae3d29cSBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
979a297a907SKarl Rupp 
9801ae3d29cSBarry Smith     *(const char**)amsopt->data = defaultv;
981a297a907SKarl Rupp 
9821ae3d29cSBarry Smith     amsopt->list  = list;
9831ae3d29cSBarry Smith     amsopt->nlist = ntext;
9841ae3d29cSBarry Smith   }
98553acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
98753acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
98853acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
98953acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
99053acd3b1SBarry Smith     }
9912aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
99253acd3b1SBarry Smith   }
99353acd3b1SBarry Smith   PetscFunctionReturn(0);
99453acd3b1SBarry Smith }
99553acd3b1SBarry Smith 
99653acd3b1SBarry Smith #undef __FUNCT__
997acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
99853acd3b1SBarry Smith /*@C
999acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1000d5649816SBarry Smith        which at most a single value can be true.
100153acd3b1SBarry Smith 
10023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100353acd3b1SBarry Smith 
100453acd3b1SBarry Smith    Input Parameters:
100553acd3b1SBarry Smith +  opt - option name
100653acd3b1SBarry Smith .  text - short string that describes the option
100753acd3b1SBarry Smith -  man - manual page with additional information on option
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Output Parameter:
101053acd3b1SBarry Smith .  flg - whether that option was set or not
101153acd3b1SBarry Smith 
101253acd3b1SBarry Smith    Level: intermediate
101353acd3b1SBarry Smith 
101453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101553acd3b1SBarry Smith 
1016acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
101753acd3b1SBarry Smith 
101853acd3b1SBarry Smith     Concepts: options database^logical group
101953acd3b1SBarry Smith 
102053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1021acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1024acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
102553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
102653acd3b1SBarry Smith @*/
10277087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
102853acd3b1SBarry Smith {
102953acd3b1SBarry Smith   PetscErrorCode ierr;
10301ae3d29cSBarry Smith   PetscOptions   amsopt;
103153acd3b1SBarry Smith 
103253acd3b1SBarry Smith   PetscFunctionBegin;
10331ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10341ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1035ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1036a297a907SKarl Rupp 
1037ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10381ae3d29cSBarry Smith   }
103968b16fdaSBarry Smith   *flg = PETSC_FALSE;
10400298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
104161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104253acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10432aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104453acd3b1SBarry Smith   }
104553acd3b1SBarry Smith   PetscFunctionReturn(0);
104653acd3b1SBarry Smith }
104753acd3b1SBarry Smith 
104853acd3b1SBarry Smith #undef __FUNCT__
1049acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
105053acd3b1SBarry Smith /*@C
1051acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1052d5649816SBarry Smith        which at most a single value can be true.
105353acd3b1SBarry Smith 
10543f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105553acd3b1SBarry Smith 
105653acd3b1SBarry Smith    Input Parameters:
105753acd3b1SBarry Smith +  opt - option name
105853acd3b1SBarry Smith .  text - short string that describes the option
105953acd3b1SBarry Smith -  man - manual page with additional information on option
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Output Parameter:
106253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith    Level: intermediate
106553acd3b1SBarry Smith 
106653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106753acd3b1SBarry Smith 
1068acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
106953acd3b1SBarry Smith 
107053acd3b1SBarry Smith     Concepts: options database^logical group
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1073acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1076acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
107753acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
107853acd3b1SBarry Smith @*/
10797087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
108053acd3b1SBarry Smith {
108153acd3b1SBarry Smith   PetscErrorCode ierr;
10821ae3d29cSBarry Smith   PetscOptions   amsopt;
108353acd3b1SBarry Smith 
108453acd3b1SBarry Smith   PetscFunctionBegin;
10851ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10861ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1087ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1088a297a907SKarl Rupp 
1089ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10901ae3d29cSBarry Smith   }
109117326d04SJed Brown   *flg = PETSC_FALSE;
10920298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10942aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109553acd3b1SBarry Smith   }
109653acd3b1SBarry Smith   PetscFunctionReturn(0);
109753acd3b1SBarry Smith }
109853acd3b1SBarry Smith 
109953acd3b1SBarry Smith #undef __FUNCT__
1100acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
110153acd3b1SBarry Smith /*@C
1102acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1103d5649816SBarry Smith        which at most a single value can be true.
110453acd3b1SBarry Smith 
11053f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110653acd3b1SBarry Smith 
110753acd3b1SBarry Smith    Input Parameters:
110853acd3b1SBarry Smith +  opt - option name
110953acd3b1SBarry Smith .  text - short string that describes the option
111053acd3b1SBarry Smith -  man - manual page with additional information on option
111153acd3b1SBarry Smith 
111253acd3b1SBarry Smith    Output Parameter:
111353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith    Level: intermediate
111653acd3b1SBarry Smith 
111753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111853acd3b1SBarry Smith 
1119acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
112053acd3b1SBarry Smith 
112153acd3b1SBarry Smith     Concepts: options database^logical group
112253acd3b1SBarry Smith 
112353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1124acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1127acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
112853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
112953acd3b1SBarry Smith @*/
11307087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
113153acd3b1SBarry Smith {
113253acd3b1SBarry Smith   PetscErrorCode ierr;
11331ae3d29cSBarry Smith   PetscOptions   amsopt;
113453acd3b1SBarry Smith 
113553acd3b1SBarry Smith   PetscFunctionBegin;
11361ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11371ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1138ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1139a297a907SKarl Rupp 
1140ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11411ae3d29cSBarry Smith   }
114217326d04SJed Brown   *flg = PETSC_FALSE;
11430298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11452aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114653acd3b1SBarry Smith   }
114753acd3b1SBarry Smith   PetscFunctionReturn(0);
114853acd3b1SBarry Smith }
114953acd3b1SBarry Smith 
115053acd3b1SBarry Smith #undef __FUNCT__
1151acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
115253acd3b1SBarry Smith /*@C
1153acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115453acd3b1SBarry Smith 
11553f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115653acd3b1SBarry Smith 
115753acd3b1SBarry Smith    Input Parameters:
115853acd3b1SBarry Smith +  opt - option name
115953acd3b1SBarry Smith .  text - short string that describes the option
116053acd3b1SBarry Smith -  man - manual page with additional information on option
116153acd3b1SBarry Smith 
116253acd3b1SBarry Smith    Output Parameter:
116353acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116453acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
116553acd3b1SBarry Smith 
116653acd3b1SBarry Smith    Level: beginner
116753acd3b1SBarry Smith 
116853acd3b1SBarry Smith    Concepts: options database^logical
116953acd3b1SBarry Smith 
117053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117153acd3b1SBarry Smith 
117253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1173acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1174acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
117553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1177acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
117853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
117953acd3b1SBarry Smith @*/
11807087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
118153acd3b1SBarry Smith {
118253acd3b1SBarry Smith   PetscErrorCode ierr;
1183ace3abfcSBarry Smith   PetscBool      iset;
1184aee2cecaSBarry Smith   PetscOptions   amsopt;
118553acd3b1SBarry Smith 
118653acd3b1SBarry Smith   PetscFunctionBegin;
1187b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1188aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1189ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1190a297a907SKarl Rupp 
1191ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1192af6d86caSBarry Smith   }
1193acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119453acd3b1SBarry Smith   if (!iset) {
119553acd3b1SBarry Smith     if (flg) *flg = deflt;
119653acd3b1SBarry Smith   }
119753acd3b1SBarry Smith   if (set) *set = iset;
119861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1199ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12002aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
120153acd3b1SBarry Smith   }
120253acd3b1SBarry Smith   PetscFunctionReturn(0);
120353acd3b1SBarry Smith }
120453acd3b1SBarry Smith 
120553acd3b1SBarry Smith #undef __FUNCT__
120653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
120753acd3b1SBarry Smith /*@C
120853acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
120953acd3b1SBarry Smith    option in the database. The values must be separated with commas with
121053acd3b1SBarry Smith    no intervening spaces.
121153acd3b1SBarry Smith 
12123f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121353acd3b1SBarry Smith 
121453acd3b1SBarry Smith    Input Parameters:
121553acd3b1SBarry Smith +  opt - the option one is seeking
121653acd3b1SBarry Smith .  text - short string describing option
121753acd3b1SBarry Smith .  man - manual page for option
121853acd3b1SBarry Smith -  nmax - maximum number of values
121953acd3b1SBarry Smith 
122053acd3b1SBarry Smith    Output Parameter:
122153acd3b1SBarry Smith +  value - location to copy values
122253acd3b1SBarry Smith .  nmax - actual number of values found
122353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122453acd3b1SBarry Smith 
122553acd3b1SBarry Smith    Level: beginner
122653acd3b1SBarry Smith 
122753acd3b1SBarry Smith    Notes:
122853acd3b1SBarry Smith    The user should pass in an array of doubles
122953acd3b1SBarry Smith 
123053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123153acd3b1SBarry Smith 
123253acd3b1SBarry Smith    Concepts: options database^array of strings
123353acd3b1SBarry Smith 
123453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1235acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1238acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
123953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
124053acd3b1SBarry Smith @*/
12417087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
124253acd3b1SBarry Smith {
124353acd3b1SBarry Smith   PetscErrorCode ierr;
124453acd3b1SBarry Smith   PetscInt       i;
1245e26ddf31SBarry Smith   PetscOptions   amsopt;
124653acd3b1SBarry Smith 
124753acd3b1SBarry Smith   PetscFunctionBegin;
1248b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1249e26ddf31SBarry Smith     PetscReal *vals;
1250e26ddf31SBarry Smith 
1251e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1252e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1253e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1254e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1255e26ddf31SBarry Smith     amsopt->arraylength = *n;
1256e26ddf31SBarry Smith   }
125753acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
125861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1259a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
126053acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1261a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
126253acd3b1SBarry Smith     }
12632aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126453acd3b1SBarry Smith   }
126553acd3b1SBarry Smith   PetscFunctionReturn(0);
126653acd3b1SBarry Smith }
126753acd3b1SBarry Smith 
126853acd3b1SBarry Smith 
126953acd3b1SBarry Smith #undef __FUNCT__
127053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
127153acd3b1SBarry Smith /*@C
127253acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1273b32a342fSShri Abhyankar    option in the database.
127453acd3b1SBarry Smith 
12753f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Input Parameters:
127853acd3b1SBarry Smith +  opt - the option one is seeking
127953acd3b1SBarry Smith .  text - short string describing option
128053acd3b1SBarry Smith .  man - manual page for option
1281f8a50e2bSBarry Smith -  n - maximum number of values
128253acd3b1SBarry Smith 
128353acd3b1SBarry Smith    Output Parameter:
128453acd3b1SBarry Smith +  value - location to copy values
1285f8a50e2bSBarry Smith .  n - actual number of values found
128653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128753acd3b1SBarry Smith 
128853acd3b1SBarry Smith    Level: beginner
128953acd3b1SBarry Smith 
129053acd3b1SBarry Smith    Notes:
1291b32a342fSShri Abhyankar    The array can be passed as
1292b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12930fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12940fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12950fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1296b32a342fSShri Abhyankar 
1297b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
129853acd3b1SBarry Smith 
129953acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
130053acd3b1SBarry Smith 
1301b32a342fSShri Abhyankar    Concepts: options database^array of ints
130253acd3b1SBarry Smith 
130353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1304acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
130553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1307acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
130853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
130953acd3b1SBarry Smith @*/
13107087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
131153acd3b1SBarry Smith {
131253acd3b1SBarry Smith   PetscErrorCode ierr;
131353acd3b1SBarry Smith   PetscInt       i;
1314e26ddf31SBarry Smith   PetscOptions   amsopt;
131553acd3b1SBarry Smith 
131653acd3b1SBarry Smith   PetscFunctionBegin;
1317b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1318e26ddf31SBarry Smith     PetscInt *vals;
1319e26ddf31SBarry Smith 
1320e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1321e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1322e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1323e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1324e26ddf31SBarry Smith     amsopt->arraylength = *n;
1325e26ddf31SBarry Smith   }
132653acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
132761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132853acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
132953acd3b1SBarry Smith     for (i=1; i<*n; i++) {
133053acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
133153acd3b1SBarry Smith     }
13322aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133353acd3b1SBarry Smith   }
133453acd3b1SBarry Smith   PetscFunctionReturn(0);
133553acd3b1SBarry Smith }
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith #undef __FUNCT__
133853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
133953acd3b1SBarry Smith /*@C
134053acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
134153acd3b1SBarry Smith    option in the database. The values must be separated with commas with
134253acd3b1SBarry Smith    no intervening spaces.
134353acd3b1SBarry Smith 
13443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
134553acd3b1SBarry Smith 
134653acd3b1SBarry Smith    Input Parameters:
134753acd3b1SBarry Smith +  opt - the option one is seeking
134853acd3b1SBarry Smith .  text - short string describing option
134953acd3b1SBarry Smith .  man - manual page for option
135053acd3b1SBarry Smith -  nmax - maximum number of strings
135153acd3b1SBarry Smith 
135253acd3b1SBarry Smith    Output Parameter:
135353acd3b1SBarry Smith +  value - location to copy strings
135453acd3b1SBarry Smith .  nmax - actual number of strings found
135553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135653acd3b1SBarry Smith 
135753acd3b1SBarry Smith    Level: beginner
135853acd3b1SBarry Smith 
135953acd3b1SBarry Smith    Notes:
136053acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
136153acd3b1SBarry Smith    strings returned by this function.
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136453acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
136553acd3b1SBarry Smith 
136653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith    Concepts: options database^array of strings
136953acd3b1SBarry Smith 
137053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1371acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1374acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
137553acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
137653acd3b1SBarry Smith @*/
13777087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
137853acd3b1SBarry Smith {
137953acd3b1SBarry Smith   PetscErrorCode ierr;
13801ae3d29cSBarry Smith   PetscOptions   amsopt;
138153acd3b1SBarry Smith 
138253acd3b1SBarry Smith   PetscFunctionBegin;
13831ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13841ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13851ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1386a297a907SKarl Rupp 
13871ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13881ae3d29cSBarry Smith   }
138953acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
139061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13912aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139253acd3b1SBarry Smith   }
139353acd3b1SBarry Smith   PetscFunctionReturn(0);
139453acd3b1SBarry Smith }
139553acd3b1SBarry Smith 
1396e2446a98SMatthew Knepley #undef __FUNCT__
1397acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1398e2446a98SMatthew Knepley /*@C
1399acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1400e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1401e2446a98SMatthew Knepley    no intervening spaces.
1402e2446a98SMatthew Knepley 
14033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1404e2446a98SMatthew Knepley 
1405e2446a98SMatthew Knepley    Input Parameters:
1406e2446a98SMatthew Knepley +  opt - the option one is seeking
1407e2446a98SMatthew Knepley .  text - short string describing option
1408e2446a98SMatthew Knepley .  man - manual page for option
1409e2446a98SMatthew Knepley -  nmax - maximum number of values
1410e2446a98SMatthew Knepley 
1411e2446a98SMatthew Knepley    Output Parameter:
1412e2446a98SMatthew Knepley +  value - location to copy values
1413e2446a98SMatthew Knepley .  nmax - actual number of values found
1414e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1415e2446a98SMatthew Knepley 
1416e2446a98SMatthew Knepley    Level: beginner
1417e2446a98SMatthew Knepley 
1418e2446a98SMatthew Knepley    Notes:
1419e2446a98SMatthew Knepley    The user should pass in an array of doubles
1420e2446a98SMatthew Knepley 
1421e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1422e2446a98SMatthew Knepley 
1423e2446a98SMatthew Knepley    Concepts: options database^array of strings
1424e2446a98SMatthew Knepley 
1425e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1426acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1427e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1428e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1429acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1430e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1431e2446a98SMatthew Knepley @*/
14327087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1433e2446a98SMatthew Knepley {
1434e2446a98SMatthew Knepley   PetscErrorCode ierr;
1435e2446a98SMatthew Knepley   PetscInt       i;
14361ae3d29cSBarry Smith   PetscOptions   amsopt;
1437e2446a98SMatthew Knepley 
1438e2446a98SMatthew Knepley   PetscFunctionBegin;
14391ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1440ace3abfcSBarry Smith     PetscBool *vals;
14411ae3d29cSBarry Smith 
14421ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL_ARRAY,&amsopt);CHKERRQ(ierr);
1443ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1444ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14451ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14461ae3d29cSBarry Smith     amsopt->arraylength = *n;
14471ae3d29cSBarry Smith   }
1448acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1449e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1450e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1451e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1452e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1453e2446a98SMatthew Knepley     }
14542aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1455e2446a98SMatthew Knepley   }
1456e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1457e2446a98SMatthew Knepley }
1458e2446a98SMatthew Knepley 
14598cc676e6SMatthew G Knepley #undef __FUNCT__
14608cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14618cc676e6SMatthew G Knepley /*@C
14628cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14638cc676e6SMatthew G Knepley 
14648cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14658cc676e6SMatthew G Knepley 
14668cc676e6SMatthew G Knepley    Input Parameters:
14678cc676e6SMatthew G Knepley +  opt - option name
14688cc676e6SMatthew G Knepley .  text - short string that describes the option
14698cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14708cc676e6SMatthew G Knepley 
14718cc676e6SMatthew G Knepley    Output Parameter:
14728cc676e6SMatthew G Knepley +  viewer - the viewer
14738cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14748cc676e6SMatthew G Knepley 
14758cc676e6SMatthew G Knepley    Level: beginner
14768cc676e6SMatthew G Knepley 
14778cc676e6SMatthew G Knepley    Concepts: options database^has int
14788cc676e6SMatthew G Knepley 
14798cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14808cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14818cc676e6SMatthew 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
14828cc676e6SMatthew G Knepley $                                     about the object to standard out
14838cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14848cc676e6SMatthew G Knepley $       draw
14858cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14868cc676e6SMatthew G Knepley 
1487cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14888cc676e6SMatthew G Knepley 
14898cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14908cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14918cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14928cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14938cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14948cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14958cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14968cc676e6SMatthew G Knepley @*/
1497cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14988cc676e6SMatthew G Knepley {
14998cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15008cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15018cc676e6SMatthew G Knepley 
15028cc676e6SMatthew G Knepley   PetscFunctionBegin;
15038cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15048cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
15058cc676e6SMatthew G Knepley     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1506a297a907SKarl Rupp 
15078cc676e6SMatthew G Knepley     *(const char**)amsopt->data = "";
15088cc676e6SMatthew G Knepley   }
1509cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15108cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15118cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15128cc676e6SMatthew G Knepley   }
15138cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15148cc676e6SMatthew G Knepley }
15158cc676e6SMatthew G Knepley 
151653acd3b1SBarry Smith 
151753acd3b1SBarry Smith #undef __FUNCT__
151853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
151953acd3b1SBarry Smith /*@C
1520b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
152153acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
152253acd3b1SBarry Smith 
15233f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
152453acd3b1SBarry Smith 
152553acd3b1SBarry Smith    Input Parameter:
152653acd3b1SBarry Smith .   head - the heading text
152753acd3b1SBarry Smith 
152853acd3b1SBarry Smith 
152953acd3b1SBarry Smith    Level: intermediate
153053acd3b1SBarry Smith 
153153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
153253acd3b1SBarry Smith 
1533b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
153453acd3b1SBarry Smith 
153553acd3b1SBarry Smith    Concepts: options database^subheading
153653acd3b1SBarry Smith 
153753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1538acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
154053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1541acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
154253acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
154353acd3b1SBarry Smith @*/
15447087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154553acd3b1SBarry Smith {
154653acd3b1SBarry Smith   PetscErrorCode ierr;
154753acd3b1SBarry Smith 
154853acd3b1SBarry Smith   PetscFunctionBegin;
154961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
155053acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
155153acd3b1SBarry Smith   }
155253acd3b1SBarry Smith   PetscFunctionReturn(0);
155353acd3b1SBarry Smith }
155453acd3b1SBarry Smith 
155553acd3b1SBarry Smith 
155653acd3b1SBarry Smith 
155753acd3b1SBarry Smith 
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith 
1560