xref: /petsc/src/sys/objects/aoptions.c (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
17d0a6c19SBarry Smith 
253acd3b1SBarry Smith /*
33fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
43fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
553acd3b1SBarry Smith 
653acd3b1SBarry Smith */
753acd3b1SBarry Smith 
8afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
9665c2dedSJed Brown #include <petscviewer.h>
1053acd3b1SBarry Smith 
112aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
122aa6d131SJed Brown 
1353acd3b1SBarry Smith /*
1453acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
153fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1653acd3b1SBarry Smith 
1753acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1853acd3b1SBarry Smith */
19f8d0b74dSMatthew Knepley PetscOptionsObjectType PetscOptionsObject;
2053acd3b1SBarry Smith PetscInt               PetscOptionsPublishCount = 0;
2153acd3b1SBarry Smith 
2253acd3b1SBarry Smith #undef __FUNCT__
2353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsBegin_Private"
2453acd3b1SBarry Smith /*
2553acd3b1SBarry Smith     Handles setting up the data structure in a call to PetscOptionsBegin()
2653acd3b1SBarry Smith */
2753acd3b1SBarry Smith PetscErrorCode PetscOptionsBegin_Private(MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
3253acd3b1SBarry Smith   PetscOptionsObject.next          = 0;
3353acd3b1SBarry Smith   PetscOptionsObject.comm          = comm;
346356e834SBarry Smith   PetscOptionsObject.changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
3753acd3b1SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject.prefix);CHKERRQ(ierr);
38c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
3953acd3b1SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject.title);CHKERRQ(ierr);
4053acd3b1SBarry Smith 
410298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject.printhelp);CHKERRQ(ierr);
4253acd3b1SBarry Smith   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1) {
4361b37b28SSatish Balay     if (!PetscOptionsObject.alreadyprinted) {
4453acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4553acd3b1SBarry Smith     }
4661b37b28SSatish Balay   }
4753acd3b1SBarry Smith   PetscFunctionReturn(0);
4853acd3b1SBarry Smith }
4953acd3b1SBarry Smith 
503194b578SJed Brown #undef __FUNCT__
513194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
523194b578SJed Brown /*
533194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
543194b578SJed Brown */
553194b578SJed Brown PetscErrorCode PetscObjectOptionsBegin_Private(PetscObject obj)
563194b578SJed Brown {
573194b578SJed Brown   PetscErrorCode ierr;
583194b578SJed Brown   char           title[256];
593194b578SJed Brown   PetscBool      flg;
603194b578SJed Brown 
613194b578SJed Brown   PetscFunctionBegin;
623194b578SJed Brown   PetscValidHeader(obj,1);
633194b578SJed Brown   PetscOptionsObject.object         = obj;
643194b578SJed Brown   PetscOptionsObject.alreadyprinted = obj->optionsprinted;
65a297a907SKarl Rupp 
663194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
673194b578SJed Brown   if (flg) {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   } else {
708caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
713194b578SJed Brown   }
723194b578SJed Brown   ierr = PetscOptionsBegin_Private(obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
733194b578SJed Brown   PetscFunctionReturn(0);
743194b578SJed Brown }
753194b578SJed Brown 
7653acd3b1SBarry Smith /*
7753acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7853acd3b1SBarry Smith */
7953acd3b1SBarry Smith #undef __FUNCT__
8053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
81e3ed6ec8SBarry Smith static int PetscOptionsCreate_Private(const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8253acd3b1SBarry Smith {
8353acd3b1SBarry Smith   int          ierr;
8453acd3b1SBarry Smith   PetscOptions next;
853be6e4c3SJed Brown   PetscBool    valid;
8653acd3b1SBarry Smith 
8753acd3b1SBarry Smith   PetscFunctionBegin;
883be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
893be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
903be6e4c3SJed Brown 
916e655a9bSJed Brown   ierr            = PetscNew(struct _n_PetscOptions,amsopt);CHKERRQ(ierr);
9253acd3b1SBarry Smith   (*amsopt)->next = 0;
9353acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
946356e834SBarry Smith   (*amsopt)->type = t;
9553acd3b1SBarry Smith   (*amsopt)->data = 0;
9661b37b28SSatish Balay 
9753acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9853acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
996356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
10053acd3b1SBarry Smith 
101a297a907SKarl Rupp   if (!PetscOptionsObject.next) PetscOptionsObject.next = *amsopt;
102a297a907SKarl Rupp   else {
10353acd3b1SBarry Smith     next = PetscOptionsObject.next;
10453acd3b1SBarry Smith     while (next->next) next = next->next;
10553acd3b1SBarry Smith     next->next = *amsopt;
10653acd3b1SBarry Smith   }
10753acd3b1SBarry Smith   PetscFunctionReturn(0);
10853acd3b1SBarry Smith }
10953acd3b1SBarry Smith 
11053acd3b1SBarry Smith #undef __FUNCT__
111aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
112aee2cecaSBarry Smith /*
1133fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith     Collective on MPI_Comm
1163fc1eb6aSBarry Smith 
1173fc1eb6aSBarry Smith    Input Parameters:
1183fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1193fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1203fc1eb6aSBarry Smith -     str - location to store input
121aee2cecaSBarry Smith 
122aee2cecaSBarry Smith     Bugs:
123aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
124aee2cecaSBarry Smith 
125aee2cecaSBarry Smith */
1263fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
127aee2cecaSBarry Smith {
128330cf3c9SBarry Smith   size_t         i;
129aee2cecaSBarry Smith   char           c;
1303fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
131aee2cecaSBarry Smith   PetscErrorCode ierr;
132aee2cecaSBarry Smith 
133aee2cecaSBarry Smith   PetscFunctionBegin;
134aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
135aee2cecaSBarry Smith   if (!rank) {
136aee2cecaSBarry Smith     c = (char) getchar();
137aee2cecaSBarry Smith     i = 0;
138aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
139aee2cecaSBarry Smith       str[i++] = c;
140aee2cecaSBarry Smith       c = (char)getchar();
141aee2cecaSBarry Smith     }
142aee2cecaSBarry Smith     str[i] = 0;
143aee2cecaSBarry Smith   }
1444dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1453fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
146aee2cecaSBarry Smith   PetscFunctionReturn(0);
147aee2cecaSBarry Smith }
148aee2cecaSBarry Smith 
149aee2cecaSBarry Smith #undef __FUNCT__
150aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
151aee2cecaSBarry Smith /*
1523cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
153aee2cecaSBarry Smith 
154aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
155aee2cecaSBarry Smith 
1567781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1577781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1587781c08eSBarry Smith 
159aee2cecaSBarry Smith     Bugs:
1607781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1613cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
162aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
163aee2cecaSBarry Smith 
1643cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1653cc1e11dSBarry Smith      address space and communicating with the PETSc program
1663cc1e11dSBarry Smith 
167aee2cecaSBarry Smith */
168aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1696356e834SBarry Smith {
1706356e834SBarry Smith   PetscErrorCode ierr;
1716356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1726356e834SBarry Smith   char           str[512];
1737781c08eSBarry Smith   PetscBool      bid;
174a4404d99SBarry Smith   PetscReal      ir,*valr;
175330cf3c9SBarry Smith   PetscInt       *vald;
176330cf3c9SBarry Smith   size_t         i;
1776356e834SBarry Smith 
178e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1796356e834SBarry Smith   while (next) {
1806356e834SBarry Smith     switch (next->type) {
1816356e834SBarry Smith     case OPTION_HEAD:
1826356e834SBarry Smith       break;
183e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
184e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
185e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
186e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
187e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
188e26ddf31SBarry Smith         if (i < next->arraylength-1) {
189e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
190e26ddf31SBarry Smith         }
191e26ddf31SBarry Smith       }
192e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
193e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
194e26ddf31SBarry Smith       if (str[0]) {
195e26ddf31SBarry Smith         PetscToken token;
196e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
197e26ddf31SBarry Smith         size_t     len;
198e26ddf31SBarry Smith         char       *value;
199ace3abfcSBarry Smith         PetscBool  foundrange;
200e26ddf31SBarry Smith 
201e26ddf31SBarry Smith         next->set = PETSC_TRUE;
202e26ddf31SBarry Smith         value     = str;
203e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
204e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
205e26ddf31SBarry Smith         while (n < nmax) {
206e26ddf31SBarry Smith           if (!value) break;
207e26ddf31SBarry Smith 
208e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
209e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
210e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
211e26ddf31SBarry Smith           if (value[0] == '-') i=2;
212e26ddf31SBarry Smith           else i=1;
213330cf3c9SBarry Smith           for (;i<len; i++) {
214e26ddf31SBarry Smith             if (value[i] == '-') {
215e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
216e26ddf31SBarry Smith               value[i] = 0;
217cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
218cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
219e32f2f54SBarry 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);
220e32f2f54SBarry 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);
221e26ddf31SBarry Smith               for (; start<end; start++) {
222e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
223e26ddf31SBarry Smith               }
224e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
225e26ddf31SBarry Smith               break;
226e26ddf31SBarry Smith             }
227e26ddf31SBarry Smith           }
228e26ddf31SBarry Smith           if (!foundrange) {
229cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
230e26ddf31SBarry Smith             dvalue++;
231e26ddf31SBarry Smith             n++;
232e26ddf31SBarry Smith           }
233e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
234e26ddf31SBarry Smith         }
2358c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
236e26ddf31SBarry Smith       }
237e26ddf31SBarry Smith       break;
238e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
239e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
240e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
241e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
242e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
243e26ddf31SBarry Smith         if (i < next->arraylength-1) {
244e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
245e26ddf31SBarry Smith         }
246e26ddf31SBarry Smith       }
247e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
248e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
249e26ddf31SBarry Smith       if (str[0]) {
250e26ddf31SBarry Smith         PetscToken token;
251e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
252e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
253e26ddf31SBarry Smith         char       *value;
254e26ddf31SBarry Smith 
255e26ddf31SBarry Smith         next->set = PETSC_TRUE;
256e26ddf31SBarry Smith         value     = str;
257e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
258e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
259e26ddf31SBarry Smith         while (n < nmax) {
260e26ddf31SBarry Smith           if (!value) break;
261cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
262e26ddf31SBarry Smith           dvalue++;
263e26ddf31SBarry Smith           n++;
264e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
265e26ddf31SBarry Smith         }
2668c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
267e26ddf31SBarry Smith       }
268e26ddf31SBarry Smith       break;
2696356e834SBarry Smith     case OPTION_INT:
270e26ddf31SBarry 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);
2713fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2723fc1eb6aSBarry Smith       if (str[0]) {
273d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
274d25d7f95SJed Brown         long long lid;
275d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
276d25d7f95SJed Brown         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %lld",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,lid);
277c272547aSJed Brown #else
278d25d7f95SJed Brown         long  lid;
279d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
280d25d7f95SJed Brown         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %ld",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,lid);
281c272547aSJed Brown #endif
282a297a907SKarl Rupp 
283d25d7f95SJed Brown         next->set = PETSC_TRUE;
284d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
285aee2cecaSBarry Smith       }
286aee2cecaSBarry Smith       break;
287aee2cecaSBarry Smith     case OPTION_REAL:
288e26ddf31SBarry 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);
2893fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2903fc1eb6aSBarry Smith       if (str[0]) {
291ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
292a4404d99SBarry Smith         sscanf(str,"%e",&ir);
293ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
294aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
295ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
296d9822059SBarry Smith         ir = strtoflt128(str,0);
297d9822059SBarry Smith #else
298513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
299a4404d99SBarry Smith #endif
300aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
301aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
302aee2cecaSBarry Smith       }
303aee2cecaSBarry Smith       break;
3047781c08eSBarry Smith     case OPTION_BOOL:
3057781c08eSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3067781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3077781c08eSBarry Smith       if (str[0]) {
3087781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3097781c08eSBarry Smith         next->set = PETSC_TRUE;
3107781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3117781c08eSBarry Smith       }
3127781c08eSBarry Smith       break;
313aee2cecaSBarry Smith     case OPTION_STRING:
314e26ddf31SBarry 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);
3153fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3163fc1eb6aSBarry Smith       if (str[0]) {
317aee2cecaSBarry Smith         next->set = PETSC_TRUE;
31877b82169SBarry Smith         ierr = PetscStrallocpy(str,(char**) &next->data);CHKERRQ(ierr);
3196356e834SBarry Smith       }
3206356e834SBarry Smith       break;
321a264d7a6SBarry Smith     case OPTION_FLIST:
322140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3233cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3243cc1e11dSBarry Smith       if (str[0]) {
3253cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3263cc1e11dSBarry Smith         next->set = PETSC_TRUE;
32777b82169SBarry Smith         ierr = PetscStrallocpy(str,(char**) &next->data);CHKERRQ(ierr);
3283cc1e11dSBarry Smith       }
3293cc1e11dSBarry Smith       break;
330b432afdaSMatthew Knepley     default:
331b432afdaSMatthew Knepley       break;
3326356e834SBarry Smith     }
3336356e834SBarry Smith     next = next->next;
3346356e834SBarry Smith   }
3356356e834SBarry Smith   PetscFunctionReturn(0);
3366356e834SBarry Smith }
3376356e834SBarry Smith 
338e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
339e04113cfSBarry Smith #include <petscviewersaws.h>
340d5649816SBarry Smith 
341d5649816SBarry Smith static int count = 0;
342d5649816SBarry Smith 
343b3506946SBarry Smith #undef __FUNCT__
344e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
345e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
346d5649816SBarry Smith {
3472657e9d9SBarry Smith   PetscFunctionBegin;
348d5649816SBarry Smith   PetscFunctionReturn(0);
349d5649816SBarry Smith }
350d5649816SBarry Smith 
351d5649816SBarry Smith #undef __FUNCT__
3527781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
353b3506946SBarry Smith /*
3547781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
355b3506946SBarry Smith 
356b3506946SBarry Smith     Bugs:
357b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
358b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
359b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
360b3506946SBarry Smith 
361b3506946SBarry Smith 
362b3506946SBarry Smith */
363475446a1SBarry Smith PetscErrorCode PetscOptionsSAWsInput()
364b3506946SBarry Smith {
365b3506946SBarry Smith   PetscErrorCode ierr;
366b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
367d5649816SBarry Smith   static int     mancount = 0;
368b3506946SBarry Smith   char           options[16];
3697aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
37088a9752cSBarry Smith   char           manname[16],textname[16];
3712657e9d9SBarry Smith   char           dir[1024];
372b3506946SBarry Smith 
373b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
374b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
375a297a907SKarl Rupp 
376e04113cfSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* SAWs will change this, so cannot pass prefix directly */
3771bc75a8dSBarry Smith 
378d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
3797781c08eSBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.title,1,SAWs_READ,SAWs_STRING));
3807781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
3812657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject.pprefix,1,SAWs_READ,SAWs_STRING));
3822657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
383b3506946SBarry Smith 
384b3506946SBarry Smith   while (next) {
385d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
3862657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
3877781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
388d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
3897781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
3907781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
3919f32e415SBarry Smith 
392b3506946SBarry Smith     switch (next->type) {
393b3506946SBarry Smith     case OPTION_HEAD:
394b3506946SBarry Smith       break;
395b3506946SBarry Smith     case OPTION_INT_ARRAY:
3967781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
3972657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
398b3506946SBarry Smith       break;
399b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4007781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4012657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
402b3506946SBarry Smith       break;
403b3506946SBarry Smith     case OPTION_INT:
4047781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4052657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
406b3506946SBarry Smith       break;
407b3506946SBarry Smith     case OPTION_REAL:
4087781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4092657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
410b3506946SBarry Smith       break;
4117781c08eSBarry Smith     case OPTION_BOOL:
4127781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4132657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4141ae3d29cSBarry Smith       break;
4157781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4167781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4172657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
41871f08665SBarry Smith       break;
419b3506946SBarry Smith     case OPTION_STRING:
4207781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4217781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4221ae3d29cSBarry Smith       break;
4231ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4247781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4252657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
426b3506946SBarry Smith       break;
427a264d7a6SBarry Smith     case OPTION_FLIST:
428a264d7a6SBarry Smith       {
429a264d7a6SBarry Smith       PetscInt ntext;
4307781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4317781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
432a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
433a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
434a264d7a6SBarry Smith       }
435a264d7a6SBarry Smith       break;
4361ae3d29cSBarry Smith     case OPTION_ELIST:
437a264d7a6SBarry Smith       {
438a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4397781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4407781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
441*785e854fSJed Brown       ierr = PetscMalloc1((ntext+1),&next->edata);CHKERRQ(ierr);
4421ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
443a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
444a264d7a6SBarry Smith       }
445a264d7a6SBarry Smith       break;
446b3506946SBarry Smith     default:
447b3506946SBarry Smith       break;
448b3506946SBarry Smith     }
449b3506946SBarry Smith     next = next->next;
450b3506946SBarry Smith   }
451b3506946SBarry Smith 
452b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4537aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
454b3506946SBarry Smith 
45588a9752cSBarry Smith   /* determine if any values have been set in GUI */
45688a9752cSBarry Smith   next = PetscOptionsObject.next;
45788a9752cSBarry Smith   while (next) {
45888a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
45988a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
46088a9752cSBarry Smith     next = next->next;
46188a9752cSBarry Smith   }
46288a9752cSBarry Smith 
463b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
464b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
465b3506946SBarry Smith 
4669a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
467b3506946SBarry Smith   PetscFunctionReturn(0);
468b3506946SBarry Smith }
469b3506946SBarry Smith #endif
470b3506946SBarry Smith 
4716356e834SBarry Smith #undef __FUNCT__
47253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
47353acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47453acd3b1SBarry Smith {
47553acd3b1SBarry Smith   PetscErrorCode ierr;
4766356e834SBarry Smith   PetscOptions   last;
4776356e834SBarry Smith   char           option[256],value[1024],tmp[32];
478330cf3c9SBarry Smith   size_t         j;
47953acd3b1SBarry Smith 
48053acd3b1SBarry Smith   PetscFunctionBegin;
481aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
482b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
483a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
484475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
485b3506946SBarry Smith #else
48671f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
487b3506946SBarry Smith #endif
488aee2cecaSBarry Smith     }
489aee2cecaSBarry Smith   }
4906356e834SBarry Smith 
491c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
492c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4936356e834SBarry Smith 
4946356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4956356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49661b37b28SSatish Balay   /* reset alreadyprinted flag */
49761b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4983194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4990298fd71SBarry Smith   PetscOptionsObject.object = NULL;
5006356e834SBarry Smith 
5016356e834SBarry Smith   while (PetscOptionsObject.next) {
5026356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5036356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5046356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5056356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5066356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5076356e834SBarry Smith       } else {
5086356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5096356e834SBarry Smith       }
5106356e834SBarry Smith 
5116356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5126356e834SBarry Smith       case OPTION_HEAD:
5136356e834SBarry Smith         break;
514e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
515e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
516e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
517e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
518e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
519e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
520e26ddf31SBarry Smith         }
521e26ddf31SBarry Smith         break;
5226356e834SBarry Smith       case OPTION_INT:
5237a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5246356e834SBarry Smith         break;
5256356e834SBarry Smith       case OPTION_REAL:
5267a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5276356e834SBarry Smith         break;
5286356e834SBarry Smith       case OPTION_REAL_ARRAY:
5297a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5306356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5317a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5326356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5336356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5346356e834SBarry Smith         }
5356356e834SBarry Smith         break;
5367781c08eSBarry Smith       case OPTION_BOOL:
53771f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5386356e834SBarry Smith         break;
5397781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
540ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5411ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
542ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5431ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5441ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5451ae3d29cSBarry Smith         }
5461ae3d29cSBarry Smith         break;
547a264d7a6SBarry Smith       case OPTION_FLIST:
5481ae3d29cSBarry Smith       case OPTION_ELIST:
5497781c08eSBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5506356e834SBarry Smith         break;
5511ae3d29cSBarry Smith       case OPTION_STRING:
552475446a1SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject.next->data);CHKERRQ(ierr);
5538da2146eSBarry Smith         break;
5541ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5551ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5561ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5571ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5581ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5591ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5601ae3d29cSBarry Smith         }
5616356e834SBarry Smith         break;
5626356e834SBarry Smith       }
5636356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5646356e834SBarry Smith     }
565503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
566503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5676356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56871f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
5697781c08eSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
5707781c08eSBarry Smith 
5716356e834SBarry Smith     last                    = PetscOptionsObject.next;
5726356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5736356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5746356e834SBarry Smith   }
5756356e834SBarry Smith   PetscOptionsObject.next = 0;
57653acd3b1SBarry Smith   PetscFunctionReturn(0);
57753acd3b1SBarry Smith }
57853acd3b1SBarry Smith 
57953acd3b1SBarry Smith #undef __FUNCT__
58053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
58153acd3b1SBarry Smith /*@C
58253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58353acd3b1SBarry Smith 
5843f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58553acd3b1SBarry Smith 
58653acd3b1SBarry Smith    Input Parameters:
58753acd3b1SBarry Smith +  opt - option name
58853acd3b1SBarry Smith .  text - short string that describes the option
58953acd3b1SBarry Smith .  man - manual page with additional information on option
59053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
59153acd3b1SBarry Smith -  defaultv - the default (current) value
59253acd3b1SBarry Smith 
59353acd3b1SBarry Smith    Output Parameter:
59453acd3b1SBarry Smith +  value - the  value to return
595b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
59653acd3b1SBarry Smith 
59753acd3b1SBarry Smith    Level: beginner
59853acd3b1SBarry Smith 
59953acd3b1SBarry Smith    Concepts: options database
60053acd3b1SBarry Smith 
60153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
60253acd3b1SBarry Smith 
60353acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60453acd3b1SBarry Smith 
60553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
606acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
607acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
610acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
611a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
61253acd3b1SBarry Smith @*/
6137087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61453acd3b1SBarry Smith {
61553acd3b1SBarry Smith   PetscErrorCode ierr;
61653acd3b1SBarry Smith   PetscInt       ntext = 0;
617aa5bb8c0SSatish Balay   PetscInt       tval;
618ace3abfcSBarry Smith   PetscBool      tflg;
61953acd3b1SBarry Smith 
62053acd3b1SBarry Smith   PetscFunctionBegin;
62153acd3b1SBarry Smith   while (list[ntext++]) {
622e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62353acd3b1SBarry Smith   }
624e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62553acd3b1SBarry Smith   ntext -= 3;
626aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
627aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
628aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
629aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
63053acd3b1SBarry Smith   PetscFunctionReturn(0);
63153acd3b1SBarry Smith }
63253acd3b1SBarry Smith 
63353acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63453acd3b1SBarry Smith #undef __FUNCT__
63553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63653acd3b1SBarry Smith /*@C
63753acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63853acd3b1SBarry Smith 
6393f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64053acd3b1SBarry Smith 
64153acd3b1SBarry Smith    Input Parameters:
64253acd3b1SBarry Smith +  opt - option name
64353acd3b1SBarry Smith .  text - short string that describes the option
64453acd3b1SBarry Smith .  man - manual page with additional information on option
64553acd3b1SBarry Smith -  defaultv - the default (current) value
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith    Output Parameter:
64853acd3b1SBarry Smith +  value - the integer value to return
64953acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Level: beginner
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Concepts: options database^has int
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
658acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
659acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
66053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
66153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
662acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
663a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
66453acd3b1SBarry Smith @*/
6657087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66653acd3b1SBarry Smith {
66753acd3b1SBarry Smith   PetscErrorCode ierr;
6686356e834SBarry Smith   PetscOptions   amsopt;
66953acd3b1SBarry Smith 
67053acd3b1SBarry Smith   PetscFunctionBegin;
671b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6726356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6736356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
674a297a907SKarl Rupp 
6756356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
676af6d86caSBarry Smith   }
67753acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6792aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
68053acd3b1SBarry Smith   }
68153acd3b1SBarry Smith   PetscFunctionReturn(0);
68253acd3b1SBarry Smith }
68353acd3b1SBarry Smith 
68453acd3b1SBarry Smith #undef __FUNCT__
68553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68653acd3b1SBarry Smith /*@C
68753acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68853acd3b1SBarry Smith 
6893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
69053acd3b1SBarry Smith 
69153acd3b1SBarry Smith    Input Parameters:
69253acd3b1SBarry Smith +  opt - option name
69353acd3b1SBarry Smith .  text - short string that describes the option
69453acd3b1SBarry Smith .  man - manual page with additional information on option
695bcbf2dc5SJed Brown .  defaultv - the default (current) value
696bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69753acd3b1SBarry Smith 
69853acd3b1SBarry Smith    Output Parameter:
69953acd3b1SBarry Smith +  value - the value to return
70053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
70153acd3b1SBarry Smith 
70253acd3b1SBarry Smith    Level: beginner
70353acd3b1SBarry Smith 
70453acd3b1SBarry Smith    Concepts: options database^has int
70553acd3b1SBarry Smith 
70653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70753acd3b1SBarry Smith 
7087fccdfe4SBarry 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).
7097fccdfe4SBarry Smith 
71053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
711acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
712acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
715acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
716a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
71753acd3b1SBarry Smith @*/
7187087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71953acd3b1SBarry Smith {
72053acd3b1SBarry Smith   PetscErrorCode ierr;
721aee2cecaSBarry Smith   PetscOptions   amsopt;
72253acd3b1SBarry Smith 
72353acd3b1SBarry Smith   PetscFunctionBegin;
724b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
725aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
72677b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
727af6d86caSBarry Smith   }
72853acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
72961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7302aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
73153acd3b1SBarry Smith   }
73253acd3b1SBarry Smith   PetscFunctionReturn(0);
73353acd3b1SBarry Smith }
73453acd3b1SBarry Smith 
73553acd3b1SBarry Smith #undef __FUNCT__
73653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsReal"
73753acd3b1SBarry Smith /*@C
73853acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
73953acd3b1SBarry Smith 
7403f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
74153acd3b1SBarry Smith 
74253acd3b1SBarry Smith    Input Parameters:
74353acd3b1SBarry Smith +  opt - option name
74453acd3b1SBarry Smith .  text - short string that describes the option
74553acd3b1SBarry Smith .  man - manual page with additional information on option
74653acd3b1SBarry Smith -  defaultv - the default (current) value
74753acd3b1SBarry Smith 
74853acd3b1SBarry Smith    Output Parameter:
74953acd3b1SBarry Smith +  value - the value to return
75053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
75153acd3b1SBarry Smith 
75253acd3b1SBarry Smith    Level: beginner
75353acd3b1SBarry Smith 
75453acd3b1SBarry Smith    Concepts: options database^has int
75553acd3b1SBarry Smith 
75653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
759acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
760acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
76153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
76253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
763acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
764a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
76553acd3b1SBarry Smith @*/
7667087cfbeSBarry Smith PetscErrorCode  PetscOptionsReal(const char opt[],const char text[],const char man[],PetscReal defaultv,PetscReal *value,PetscBool  *set)
76753acd3b1SBarry Smith {
76853acd3b1SBarry Smith   PetscErrorCode ierr;
769538aa990SBarry Smith   PetscOptions   amsopt;
77053acd3b1SBarry Smith 
77153acd3b1SBarry Smith   PetscFunctionBegin;
772b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
773538aa990SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
774538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
775a297a907SKarl Rupp 
776538aa990SBarry Smith     *(PetscReal*)amsopt->data = defaultv;
777538aa990SBarry Smith   }
77853acd3b1SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
77961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
7802aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
78153acd3b1SBarry Smith   }
78253acd3b1SBarry Smith   PetscFunctionReturn(0);
78353acd3b1SBarry Smith }
78453acd3b1SBarry Smith 
78553acd3b1SBarry Smith #undef __FUNCT__
78653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsScalar"
78753acd3b1SBarry Smith /*@C
78853acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
78953acd3b1SBarry Smith 
7903f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
79153acd3b1SBarry Smith 
79253acd3b1SBarry Smith    Input Parameters:
79353acd3b1SBarry Smith +  opt - option name
79453acd3b1SBarry Smith .  text - short string that describes the option
79553acd3b1SBarry Smith .  man - manual page with additional information on option
79653acd3b1SBarry Smith -  defaultv - the default (current) value
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith    Output Parameter:
79953acd3b1SBarry Smith +  value - the value to return
80053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
80153acd3b1SBarry Smith 
80253acd3b1SBarry Smith    Level: beginner
80353acd3b1SBarry Smith 
80453acd3b1SBarry Smith    Concepts: options database^has int
80553acd3b1SBarry Smith 
80653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
80753acd3b1SBarry Smith 
80853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
809acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
810acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
81153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
81253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
813acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
814a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
81553acd3b1SBarry Smith @*/
8167087cfbeSBarry Smith PetscErrorCode  PetscOptionsScalar(const char opt[],const char text[],const char man[],PetscScalar defaultv,PetscScalar *value,PetscBool  *set)
81753acd3b1SBarry Smith {
81853acd3b1SBarry Smith   PetscErrorCode ierr;
81953acd3b1SBarry Smith 
82053acd3b1SBarry Smith   PetscFunctionBegin;
82153acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
82253acd3b1SBarry Smith   ierr = PetscOptionsReal(opt,text,man,defaultv,value,set);CHKERRQ(ierr);
82353acd3b1SBarry Smith #else
82453acd3b1SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
82553acd3b1SBarry Smith #endif
82653acd3b1SBarry Smith   PetscFunctionReturn(0);
82753acd3b1SBarry Smith }
82853acd3b1SBarry Smith 
82953acd3b1SBarry Smith #undef __FUNCT__
83053acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsName"
83153acd3b1SBarry Smith /*@C
83290d69ab7SBarry 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
83390d69ab7SBarry Smith                       its value is set to false.
83453acd3b1SBarry Smith 
8353f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83653acd3b1SBarry Smith 
83753acd3b1SBarry Smith    Input Parameters:
83853acd3b1SBarry Smith +  opt - option name
83953acd3b1SBarry Smith .  text - short string that describes the option
84053acd3b1SBarry Smith -  man - manual page with additional information on option
84153acd3b1SBarry Smith 
84253acd3b1SBarry Smith    Output Parameter:
84353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
84453acd3b1SBarry Smith 
84553acd3b1SBarry Smith    Level: beginner
84653acd3b1SBarry Smith 
84753acd3b1SBarry Smith    Concepts: options database^has int
84853acd3b1SBarry Smith 
84953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
852acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
853acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
856acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
857a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
85853acd3b1SBarry Smith @*/
8597087cfbeSBarry Smith PetscErrorCode  PetscOptionsName(const char opt[],const char text[],const char man[],PetscBool  *flg)
86053acd3b1SBarry Smith {
86153acd3b1SBarry Smith   PetscErrorCode ierr;
8621ae3d29cSBarry Smith   PetscOptions   amsopt;
86353acd3b1SBarry Smith 
86453acd3b1SBarry Smith   PetscFunctionBegin;
8651ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
8667781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
867ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
868a297a907SKarl Rupp 
869ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
8701ae3d29cSBarry Smith   }
87153acd3b1SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject.prefix,opt,flg);CHKERRQ(ierr);
87261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
8732aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
87453acd3b1SBarry Smith   }
87553acd3b1SBarry Smith   PetscFunctionReturn(0);
87653acd3b1SBarry Smith }
87753acd3b1SBarry Smith 
87853acd3b1SBarry Smith #undef __FUNCT__
879a264d7a6SBarry Smith #define __FUNCT__ "PetscOptionsFList"
88053acd3b1SBarry Smith /*@C
881a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
88253acd3b1SBarry Smith 
8833f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88453acd3b1SBarry Smith 
88553acd3b1SBarry Smith    Input Parameters:
88653acd3b1SBarry Smith +  opt - option name
88753acd3b1SBarry Smith .  text - short string that describes the option
88853acd3b1SBarry Smith .  man - manual page with additional information on option
88953acd3b1SBarry Smith .  list - the possible choices
8903cc1e11dSBarry Smith .  defaultv - the default (current) value
8913cc1e11dSBarry Smith -  len - the length of the character array value
89253acd3b1SBarry Smith 
89353acd3b1SBarry Smith    Output Parameter:
89453acd3b1SBarry Smith +  value - the value to return
89553acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
89653acd3b1SBarry Smith 
89753acd3b1SBarry Smith    Level: intermediate
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
90053acd3b1SBarry Smith 
90153acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
90253acd3b1SBarry Smith 
90353acd3b1SBarry Smith    To get a listing of all currently specified options,
90488c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
90553acd3b1SBarry Smith 
906eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
907eabe10d7SBarry Smith 
90853acd3b1SBarry Smith    Concepts: options database^list
90953acd3b1SBarry Smith 
91053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
911acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
91253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
914acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
915a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
91653acd3b1SBarry Smith @*/
917a264d7a6SBarry Smith PetscErrorCode  PetscOptionsFList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91853acd3b1SBarry Smith {
91953acd3b1SBarry Smith   PetscErrorCode ierr;
9203cc1e11dSBarry Smith   PetscOptions   amsopt;
92153acd3b1SBarry Smith 
92253acd3b1SBarry Smith   PetscFunctionBegin;
923b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
924a264d7a6SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
92577b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
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
946a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
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 
958a264d7a6SBarry Smith    See PetscOptionsFList() 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(),
967a264d7a6SBarry Smith           PetscOptionsFList(), 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);
97877b82169SBarry Smith     ierr = PetscStrallocpy(defaultv ? defaultv : "",(char**) &amsopt->data);CHKERRQ(ierr);
9791ae3d29cSBarry Smith     amsopt->list  = list;
9801ae3d29cSBarry Smith     amsopt->nlist = ntext;
9811ae3d29cSBarry Smith   }
98253acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98361b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
98453acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
98553acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
98653acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
98753acd3b1SBarry Smith     }
9882aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
98953acd3b1SBarry Smith   }
99053acd3b1SBarry Smith   PetscFunctionReturn(0);
99153acd3b1SBarry Smith }
99253acd3b1SBarry Smith 
99353acd3b1SBarry Smith #undef __FUNCT__
994acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
99553acd3b1SBarry Smith /*@C
996acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
997d5649816SBarry Smith        which at most a single value can be true.
99853acd3b1SBarry Smith 
9993f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100053acd3b1SBarry Smith 
100153acd3b1SBarry Smith    Input Parameters:
100253acd3b1SBarry Smith +  opt - option name
100353acd3b1SBarry Smith .  text - short string that describes the option
100453acd3b1SBarry Smith -  man - manual page with additional information on option
100553acd3b1SBarry Smith 
100653acd3b1SBarry Smith    Output Parameter:
100753acd3b1SBarry Smith .  flg - whether that option was set or not
100853acd3b1SBarry Smith 
100953acd3b1SBarry Smith    Level: intermediate
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101253acd3b1SBarry Smith 
1013acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith     Concepts: options database^logical group
101653acd3b1SBarry Smith 
101753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1018acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1021acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1022a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
102353acd3b1SBarry Smith @*/
10247087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
102553acd3b1SBarry Smith {
102653acd3b1SBarry Smith   PetscErrorCode ierr;
10271ae3d29cSBarry Smith   PetscOptions   amsopt;
102853acd3b1SBarry Smith 
102953acd3b1SBarry Smith   PetscFunctionBegin;
10301ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10317781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1032ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1033a297a907SKarl Rupp 
1034ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10351ae3d29cSBarry Smith   }
103668b16fdaSBarry Smith   *flg = PETSC_FALSE;
10370298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
103861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
103953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10402aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104153acd3b1SBarry Smith   }
104253acd3b1SBarry Smith   PetscFunctionReturn(0);
104353acd3b1SBarry Smith }
104453acd3b1SBarry Smith 
104553acd3b1SBarry Smith #undef __FUNCT__
1046acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
104753acd3b1SBarry Smith /*@C
1048acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1049d5649816SBarry Smith        which at most a single value can be true.
105053acd3b1SBarry Smith 
10513f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105253acd3b1SBarry Smith 
105353acd3b1SBarry Smith    Input Parameters:
105453acd3b1SBarry Smith +  opt - option name
105553acd3b1SBarry Smith .  text - short string that describes the option
105653acd3b1SBarry Smith -  man - manual page with additional information on option
105753acd3b1SBarry Smith 
105853acd3b1SBarry Smith    Output Parameter:
105953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Level: intermediate
106253acd3b1SBarry Smith 
106353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106453acd3b1SBarry Smith 
1065acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith     Concepts: options database^logical group
106853acd3b1SBarry Smith 
106953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1070acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1073acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1074a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107553acd3b1SBarry Smith @*/
10767087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
107753acd3b1SBarry Smith {
107853acd3b1SBarry Smith   PetscErrorCode ierr;
10791ae3d29cSBarry Smith   PetscOptions   amsopt;
108053acd3b1SBarry Smith 
108153acd3b1SBarry Smith   PetscFunctionBegin;
10821ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10837781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1084ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1085a297a907SKarl Rupp 
1086ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10871ae3d29cSBarry Smith   }
108817326d04SJed Brown   *flg = PETSC_FALSE;
10890298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10912aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109253acd3b1SBarry Smith   }
109353acd3b1SBarry Smith   PetscFunctionReturn(0);
109453acd3b1SBarry Smith }
109553acd3b1SBarry Smith 
109653acd3b1SBarry Smith #undef __FUNCT__
1097acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
109853acd3b1SBarry Smith /*@C
1099acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1100d5649816SBarry Smith        which at most a single value can be true.
110153acd3b1SBarry Smith 
11023f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110353acd3b1SBarry Smith 
110453acd3b1SBarry Smith    Input Parameters:
110553acd3b1SBarry Smith +  opt - option name
110653acd3b1SBarry Smith .  text - short string that describes the option
110753acd3b1SBarry Smith -  man - manual page with additional information on option
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith    Output Parameter:
111053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111153acd3b1SBarry Smith 
111253acd3b1SBarry Smith    Level: intermediate
111353acd3b1SBarry Smith 
111453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111553acd3b1SBarry Smith 
1116acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
111753acd3b1SBarry Smith 
111853acd3b1SBarry Smith     Concepts: options database^logical group
111953acd3b1SBarry Smith 
112053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1121acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1124acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1125a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
112653acd3b1SBarry Smith @*/
11277087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
112853acd3b1SBarry Smith {
112953acd3b1SBarry Smith   PetscErrorCode ierr;
11301ae3d29cSBarry Smith   PetscOptions   amsopt;
113153acd3b1SBarry Smith 
113253acd3b1SBarry Smith   PetscFunctionBegin;
11331ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11347781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1135ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1136a297a907SKarl Rupp 
1137ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11381ae3d29cSBarry Smith   }
113917326d04SJed Brown   *flg = PETSC_FALSE;
11400298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11422aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114353acd3b1SBarry Smith   }
114453acd3b1SBarry Smith   PetscFunctionReturn(0);
114553acd3b1SBarry Smith }
114653acd3b1SBarry Smith 
114753acd3b1SBarry Smith #undef __FUNCT__
1148acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
114953acd3b1SBarry Smith /*@C
1150acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115153acd3b1SBarry Smith 
11523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115353acd3b1SBarry Smith 
115453acd3b1SBarry Smith    Input Parameters:
115553acd3b1SBarry Smith +  opt - option name
115653acd3b1SBarry Smith .  text - short string that describes the option
115753acd3b1SBarry Smith -  man - manual page with additional information on option
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith    Output Parameter:
116053acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116153acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
116253acd3b1SBarry Smith 
116353acd3b1SBarry Smith    Level: beginner
116453acd3b1SBarry Smith 
116553acd3b1SBarry Smith    Concepts: options database^logical
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116853acd3b1SBarry Smith 
116953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1170acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1171acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
117253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1174acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1175a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
117653acd3b1SBarry Smith @*/
11777087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
117853acd3b1SBarry Smith {
117953acd3b1SBarry Smith   PetscErrorCode ierr;
1180ace3abfcSBarry Smith   PetscBool      iset;
1181aee2cecaSBarry Smith   PetscOptions   amsopt;
118253acd3b1SBarry Smith 
118353acd3b1SBarry Smith   PetscFunctionBegin;
1184b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
11857781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1186ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1187a297a907SKarl Rupp 
1188ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1189af6d86caSBarry Smith   }
1190acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119153acd3b1SBarry Smith   if (!iset) {
119253acd3b1SBarry Smith     if (flg) *flg = deflt;
119353acd3b1SBarry Smith   }
119453acd3b1SBarry Smith   if (set) *set = iset;
119561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1196ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
11972aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
119853acd3b1SBarry Smith   }
119953acd3b1SBarry Smith   PetscFunctionReturn(0);
120053acd3b1SBarry Smith }
120153acd3b1SBarry Smith 
120253acd3b1SBarry Smith #undef __FUNCT__
120353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
120453acd3b1SBarry Smith /*@C
120553acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
120653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
120753acd3b1SBarry Smith    no intervening spaces.
120853acd3b1SBarry Smith 
12093f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121053acd3b1SBarry Smith 
121153acd3b1SBarry Smith    Input Parameters:
121253acd3b1SBarry Smith +  opt - the option one is seeking
121353acd3b1SBarry Smith .  text - short string describing option
121453acd3b1SBarry Smith .  man - manual page for option
121553acd3b1SBarry Smith -  nmax - maximum number of values
121653acd3b1SBarry Smith 
121753acd3b1SBarry Smith    Output Parameter:
121853acd3b1SBarry Smith +  value - location to copy values
121953acd3b1SBarry Smith .  nmax - actual number of values found
122053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122153acd3b1SBarry Smith 
122253acd3b1SBarry Smith    Level: beginner
122353acd3b1SBarry Smith 
122453acd3b1SBarry Smith    Notes:
122553acd3b1SBarry Smith    The user should pass in an array of doubles
122653acd3b1SBarry Smith 
122753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122853acd3b1SBarry Smith 
122953acd3b1SBarry Smith    Concepts: options database^array of strings
123053acd3b1SBarry Smith 
123153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1232acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1235acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1236a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
123753acd3b1SBarry Smith @*/
12387087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
123953acd3b1SBarry Smith {
124053acd3b1SBarry Smith   PetscErrorCode ierr;
124153acd3b1SBarry Smith   PetscInt       i;
1242e26ddf31SBarry Smith   PetscOptions   amsopt;
124353acd3b1SBarry Smith 
124453acd3b1SBarry Smith   PetscFunctionBegin;
1245b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1246e26ddf31SBarry Smith     PetscReal *vals;
1247e26ddf31SBarry Smith 
1248e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1249*785e854fSJed Brown     ierr = PetscMalloc1((*n),&amsopt->data);CHKERRQ(ierr);
1250e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1251e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1252e26ddf31SBarry Smith     amsopt->arraylength = *n;
1253e26ddf31SBarry Smith   }
125453acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
125561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1256a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
125753acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1258a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
125953acd3b1SBarry Smith     }
12602aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126153acd3b1SBarry Smith   }
126253acd3b1SBarry Smith   PetscFunctionReturn(0);
126353acd3b1SBarry Smith }
126453acd3b1SBarry Smith 
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith #undef __FUNCT__
126753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
126853acd3b1SBarry Smith /*@C
126953acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1270b32a342fSShri Abhyankar    option in the database.
127153acd3b1SBarry Smith 
12723f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Input Parameters:
127553acd3b1SBarry Smith +  opt - the option one is seeking
127653acd3b1SBarry Smith .  text - short string describing option
127753acd3b1SBarry Smith .  man - manual page for option
1278f8a50e2bSBarry Smith -  n - maximum number of values
127953acd3b1SBarry Smith 
128053acd3b1SBarry Smith    Output Parameter:
128153acd3b1SBarry Smith +  value - location to copy values
1282f8a50e2bSBarry Smith .  n - actual number of values found
128353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128453acd3b1SBarry Smith 
128553acd3b1SBarry Smith    Level: beginner
128653acd3b1SBarry Smith 
128753acd3b1SBarry Smith    Notes:
1288b32a342fSShri Abhyankar    The array can be passed as
1289b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12900fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12910fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12920fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1293b32a342fSShri Abhyankar 
1294b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
129553acd3b1SBarry Smith 
129653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
129753acd3b1SBarry Smith 
1298b32a342fSShri Abhyankar    Concepts: options database^array of ints
129953acd3b1SBarry Smith 
130053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1301acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
130253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1304acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1305a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
130653acd3b1SBarry Smith @*/
13077087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
130853acd3b1SBarry Smith {
130953acd3b1SBarry Smith   PetscErrorCode ierr;
131053acd3b1SBarry Smith   PetscInt       i;
1311e26ddf31SBarry Smith   PetscOptions   amsopt;
131253acd3b1SBarry Smith 
131353acd3b1SBarry Smith   PetscFunctionBegin;
1314b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1315e26ddf31SBarry Smith     PetscInt *vals;
1316e26ddf31SBarry Smith 
1317e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1318*785e854fSJed Brown     ierr = PetscMalloc1((*n),&amsopt->data);CHKERRQ(ierr);
1319e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1320e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1321e26ddf31SBarry Smith     amsopt->arraylength = *n;
1322e26ddf31SBarry Smith   }
132353acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
132461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
132653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
132753acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
132853acd3b1SBarry Smith     }
13292aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133053acd3b1SBarry Smith   }
133153acd3b1SBarry Smith   PetscFunctionReturn(0);
133253acd3b1SBarry Smith }
133353acd3b1SBarry Smith 
133453acd3b1SBarry Smith #undef __FUNCT__
133553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
133653acd3b1SBarry Smith /*@C
133753acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
133853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
133953acd3b1SBarry Smith    no intervening spaces.
134053acd3b1SBarry Smith 
13413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
134253acd3b1SBarry Smith 
134353acd3b1SBarry Smith    Input Parameters:
134453acd3b1SBarry Smith +  opt - the option one is seeking
134553acd3b1SBarry Smith .  text - short string describing option
134653acd3b1SBarry Smith .  man - manual page for option
134753acd3b1SBarry Smith -  nmax - maximum number of strings
134853acd3b1SBarry Smith 
134953acd3b1SBarry Smith    Output Parameter:
135053acd3b1SBarry Smith +  value - location to copy strings
135153acd3b1SBarry Smith .  nmax - actual number of strings found
135253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135353acd3b1SBarry Smith 
135453acd3b1SBarry Smith    Level: beginner
135553acd3b1SBarry Smith 
135653acd3b1SBarry Smith    Notes:
135753acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
135853acd3b1SBarry Smith    strings returned by this function.
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136153acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136453acd3b1SBarry Smith 
136553acd3b1SBarry Smith    Concepts: options database^array of strings
136653acd3b1SBarry Smith 
136753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1368acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
136953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1371acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1372a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
137353acd3b1SBarry Smith @*/
13747087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
137553acd3b1SBarry Smith {
137653acd3b1SBarry Smith   PetscErrorCode ierr;
13771ae3d29cSBarry Smith   PetscOptions   amsopt;
137853acd3b1SBarry Smith 
137953acd3b1SBarry Smith   PetscFunctionBegin;
13801ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13811ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1382*785e854fSJed Brown     ierr = PetscMalloc1((*nmax),&amsopt->data);CHKERRQ(ierr);
1383a297a907SKarl Rupp 
13841ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13851ae3d29cSBarry Smith   }
138653acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
138761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13882aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
138953acd3b1SBarry Smith   }
139053acd3b1SBarry Smith   PetscFunctionReturn(0);
139153acd3b1SBarry Smith }
139253acd3b1SBarry Smith 
1393e2446a98SMatthew Knepley #undef __FUNCT__
1394acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1395e2446a98SMatthew Knepley /*@C
1396acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1397e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1398e2446a98SMatthew Knepley    no intervening spaces.
1399e2446a98SMatthew Knepley 
14003f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1401e2446a98SMatthew Knepley 
1402e2446a98SMatthew Knepley    Input Parameters:
1403e2446a98SMatthew Knepley +  opt - the option one is seeking
1404e2446a98SMatthew Knepley .  text - short string describing option
1405e2446a98SMatthew Knepley .  man - manual page for option
1406e2446a98SMatthew Knepley -  nmax - maximum number of values
1407e2446a98SMatthew Knepley 
1408e2446a98SMatthew Knepley    Output Parameter:
1409e2446a98SMatthew Knepley +  value - location to copy values
1410e2446a98SMatthew Knepley .  nmax - actual number of values found
1411e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1412e2446a98SMatthew Knepley 
1413e2446a98SMatthew Knepley    Level: beginner
1414e2446a98SMatthew Knepley 
1415e2446a98SMatthew Knepley    Notes:
1416e2446a98SMatthew Knepley    The user should pass in an array of doubles
1417e2446a98SMatthew Knepley 
1418e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1419e2446a98SMatthew Knepley 
1420e2446a98SMatthew Knepley    Concepts: options database^array of strings
1421e2446a98SMatthew Knepley 
1422e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1423acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1424e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1425e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1426acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1427a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1428e2446a98SMatthew Knepley @*/
14297087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1430e2446a98SMatthew Knepley {
1431e2446a98SMatthew Knepley   PetscErrorCode ierr;
1432e2446a98SMatthew Knepley   PetscInt       i;
14331ae3d29cSBarry Smith   PetscOptions   amsopt;
1434e2446a98SMatthew Knepley 
1435e2446a98SMatthew Knepley   PetscFunctionBegin;
14361ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1437ace3abfcSBarry Smith     PetscBool *vals;
14381ae3d29cSBarry Smith 
14397781c08eSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1440*785e854fSJed Brown     ierr = PetscMalloc1((*n),&amsopt->data);CHKERRQ(ierr);
1441ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14421ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14431ae3d29cSBarry Smith     amsopt->arraylength = *n;
14441ae3d29cSBarry Smith   }
1445acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1446e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1447e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1448e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1449e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1450e2446a98SMatthew Knepley     }
14512aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1452e2446a98SMatthew Knepley   }
1453e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1454e2446a98SMatthew Knepley }
1455e2446a98SMatthew Knepley 
14568cc676e6SMatthew G Knepley #undef __FUNCT__
14578cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14588cc676e6SMatthew G Knepley /*@C
14598cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14608cc676e6SMatthew G Knepley 
14618cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14628cc676e6SMatthew G Knepley 
14638cc676e6SMatthew G Knepley    Input Parameters:
14648cc676e6SMatthew G Knepley +  opt - option name
14658cc676e6SMatthew G Knepley .  text - short string that describes the option
14668cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14678cc676e6SMatthew G Knepley 
14688cc676e6SMatthew G Knepley    Output Parameter:
14698cc676e6SMatthew G Knepley +  viewer - the viewer
14708cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14718cc676e6SMatthew G Knepley 
14728cc676e6SMatthew G Knepley    Level: beginner
14738cc676e6SMatthew G Knepley 
14748cc676e6SMatthew G Knepley    Concepts: options database^has int
14758cc676e6SMatthew G Knepley 
14768cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14778cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14788cc676e6SMatthew 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
14798cc676e6SMatthew G Knepley $                                     about the object to standard out
14808cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14818cc676e6SMatthew G Knepley $       draw
14828cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14838cc676e6SMatthew G Knepley 
1484cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14858cc676e6SMatthew G Knepley 
14868cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14878cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14888cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14898cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14908cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14918cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1492a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
14938cc676e6SMatthew G Knepley @*/
1494cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14958cc676e6SMatthew G Knepley {
14968cc676e6SMatthew G Knepley   PetscErrorCode ierr;
14978cc676e6SMatthew G Knepley   PetscOptions   amsopt;
14988cc676e6SMatthew G Knepley 
14998cc676e6SMatthew G Knepley   PetscFunctionBegin;
15008cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15018cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
150277b82169SBarry Smith     ierr = PetscStrallocpy("",(char**) &amsopt->data);CHKERRQ(ierr);
15038cc676e6SMatthew G Knepley   }
1504cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15058cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15068cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15078cc676e6SMatthew G Knepley   }
15088cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15098cc676e6SMatthew G Knepley }
15108cc676e6SMatthew G Knepley 
151153acd3b1SBarry Smith 
151253acd3b1SBarry Smith #undef __FUNCT__
151353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
151453acd3b1SBarry Smith /*@C
1515b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
151653acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
151753acd3b1SBarry Smith 
15183f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
151953acd3b1SBarry Smith 
152053acd3b1SBarry Smith    Input Parameter:
152153acd3b1SBarry Smith .   head - the heading text
152253acd3b1SBarry Smith 
152353acd3b1SBarry Smith 
152453acd3b1SBarry Smith    Level: intermediate
152553acd3b1SBarry Smith 
152653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152753acd3b1SBarry Smith 
1528b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
152953acd3b1SBarry Smith 
153053acd3b1SBarry Smith    Concepts: options database^subheading
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1533acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
153453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
153553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1536acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1537a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
153853acd3b1SBarry Smith @*/
15397087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154053acd3b1SBarry Smith {
154153acd3b1SBarry Smith   PetscErrorCode ierr;
154253acd3b1SBarry Smith 
154353acd3b1SBarry Smith   PetscFunctionBegin;
154461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
154553acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
154653acd3b1SBarry Smith   }
154753acd3b1SBarry Smith   PetscFunctionReturn(0);
154853acd3b1SBarry Smith }
154953acd3b1SBarry Smith 
155053acd3b1SBarry Smith 
155153acd3b1SBarry Smith 
155253acd3b1SBarry Smith 
155353acd3b1SBarry Smith 
155453acd3b1SBarry Smith 
1555