xref: /petsc/src/sys/objects/aoptions.c (revision 665c2ded495bb9782a7454dcfef3abf1536c3670)
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*/
9*665c2dedSJed 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 
156aee2cecaSBarry Smith     Bugs:
157aee2cecaSBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
1583cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
159aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
160aee2cecaSBarry Smith 
1613cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1623cc1e11dSBarry Smith      address space and communicating with the PETSc program
1633cc1e11dSBarry Smith 
164aee2cecaSBarry Smith */
165aee2cecaSBarry Smith PetscErrorCode PetscOptionsGetFromTextInput()
1666356e834SBarry Smith {
1676356e834SBarry Smith   PetscErrorCode ierr;
1686356e834SBarry Smith   PetscOptions   next = PetscOptionsObject.next;
1696356e834SBarry Smith   char           str[512];
170a4404d99SBarry Smith   PetscInt       id;
171a4404d99SBarry Smith   PetscReal      ir,*valr;
172330cf3c9SBarry Smith   PetscInt       *vald;
173330cf3c9SBarry Smith   size_t         i;
1746356e834SBarry Smith 
175e26ddf31SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject.title);CHKERRQ(ierr);
1766356e834SBarry Smith   while (next) {
1776356e834SBarry Smith     switch (next->type) {
1786356e834SBarry Smith     case OPTION_HEAD:
1796356e834SBarry Smith       break;
180e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
181e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
182e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
183e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
184e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
185e26ddf31SBarry Smith         if (i < next->arraylength-1) {
186e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
187e26ddf31SBarry Smith         }
188e26ddf31SBarry Smith       }
189e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
190e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
191e26ddf31SBarry Smith       if (str[0]) {
192e26ddf31SBarry Smith         PetscToken token;
193e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
194e26ddf31SBarry Smith         size_t     len;
195e26ddf31SBarry Smith         char       *value;
196ace3abfcSBarry Smith         PetscBool  foundrange;
197e26ddf31SBarry Smith 
198e26ddf31SBarry Smith         next->set = PETSC_TRUE;
199e26ddf31SBarry Smith         value     = str;
200e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
201e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
202e26ddf31SBarry Smith         while (n < nmax) {
203e26ddf31SBarry Smith           if (!value) break;
204e26ddf31SBarry Smith 
205e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
206e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
207e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
208e26ddf31SBarry Smith           if (value[0] == '-') i=2;
209e26ddf31SBarry Smith           else i=1;
210330cf3c9SBarry Smith           for (;i<len; i++) {
211e26ddf31SBarry Smith             if (value[i] == '-') {
212e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
213e26ddf31SBarry Smith               value[i] = 0;
214cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
215cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
216e32f2f54SBarry 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);
217e32f2f54SBarry 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);
218e26ddf31SBarry Smith               for (; start<end; start++) {
219e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
220e26ddf31SBarry Smith               }
221e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
222e26ddf31SBarry Smith               break;
223e26ddf31SBarry Smith             }
224e26ddf31SBarry Smith           }
225e26ddf31SBarry Smith           if (!foundrange) {
226cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
227e26ddf31SBarry Smith             dvalue++;
228e26ddf31SBarry Smith             n++;
229e26ddf31SBarry Smith           }
230e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
231e26ddf31SBarry Smith         }
2328c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
233e26ddf31SBarry Smith       }
234e26ddf31SBarry Smith       break;
235e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
236e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",next->option+1);CHKERRQ(ierr);
237e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
238e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
239e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
240e26ddf31SBarry Smith         if (i < next->arraylength-1) {
241e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
242e26ddf31SBarry Smith         }
243e26ddf31SBarry Smith       }
244e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s)",next->text,next->man);CHKERRQ(ierr);
245e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
246e26ddf31SBarry Smith       if (str[0]) {
247e26ddf31SBarry Smith         PetscToken token;
248e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
249e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
250e26ddf31SBarry Smith         char       *value;
251e26ddf31SBarry Smith 
252e26ddf31SBarry Smith         next->set = PETSC_TRUE;
253e26ddf31SBarry Smith         value     = str;
254e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
255e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
256e26ddf31SBarry Smith         while (n < nmax) {
257e26ddf31SBarry Smith           if (!value) break;
258cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
259e26ddf31SBarry Smith           dvalue++;
260e26ddf31SBarry Smith           n++;
261e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
262e26ddf31SBarry Smith         }
2638c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
264e26ddf31SBarry Smith       }
265e26ddf31SBarry Smith       break;
2666356e834SBarry Smith     case OPTION_INT:
267e26ddf31SBarry 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);
2683fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2693fc1eb6aSBarry Smith       if (str[0]) {
270c272547aSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
271c272547aSJed Brown         sscanf(str,"%lld",&id);
272c272547aSJed Brown #else
273aee2cecaSBarry Smith         sscanf(str,"%d",&id);
274c272547aSJed Brown #endif
275aee2cecaSBarry Smith         next->set = PETSC_TRUE;
276a297a907SKarl Rupp 
277aee2cecaSBarry Smith         *((PetscInt*)next->data) = id;
278aee2cecaSBarry Smith       }
279aee2cecaSBarry Smith       break;
280aee2cecaSBarry Smith     case OPTION_REAL:
281e26ddf31SBarry 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);
2823fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2833fc1eb6aSBarry Smith       if (str[0]) {
284ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
285a4404d99SBarry Smith         sscanf(str,"%e",&ir);
286ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
287aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
288ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
289d9822059SBarry Smith         ir = strtoflt128(str,0);
290d9822059SBarry Smith #else
291513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
292a4404d99SBarry Smith #endif
293aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
294aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
295aee2cecaSBarry Smith       }
296aee2cecaSBarry Smith       break;
297aee2cecaSBarry Smith     case OPTION_LOGICAL:
298aee2cecaSBarry Smith     case OPTION_STRING:
299e26ddf31SBarry 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);
3003fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3013fc1eb6aSBarry Smith       if (str[0]) {
302aee2cecaSBarry Smith         next->set = PETSC_TRUE;
303a297a907SKarl Rupp 
304330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3056356e834SBarry Smith       }
3066356e834SBarry Smith       break;
3073cc1e11dSBarry Smith     case OPTION_LIST:
308140e18c1SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject.prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3093cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3103cc1e11dSBarry Smith       if (str[0]) {
3113cc1e11dSBarry Smith         PetscOptionsObject.changedmethod = PETSC_TRUE;
3123cc1e11dSBarry Smith         next->set = PETSC_TRUE;
313330cf3c9SBarry Smith         ierr = PetscStrcpy((char*)next->data,str);CHKERRQ(ierr);
3143cc1e11dSBarry Smith       }
3153cc1e11dSBarry Smith       break;
316b432afdaSMatthew Knepley     default:
317b432afdaSMatthew Knepley       break;
3186356e834SBarry Smith     }
3196356e834SBarry Smith     next = next->next;
3206356e834SBarry Smith   }
3216356e834SBarry Smith   PetscFunctionReturn(0);
3226356e834SBarry Smith }
3236356e834SBarry Smith 
324b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
3251ae3d29cSBarry Smith #define CHKERRAMS(err)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"AMS Error: %s",msg);}
3261ae3d29cSBarry Smith #define CHKERRAMSFieldName(err,fn)  if (err) {char *msg; AMS_Explain_error((err), &(msg)); SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Fieldname %s, AMS Error: %s",fn,msg);}
327d5649816SBarry Smith 
328d5649816SBarry Smith static int count = 0;
329d5649816SBarry Smith 
330b3506946SBarry Smith #undef __FUNCT__
331d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSDestroy"
332d5649816SBarry Smith PetscErrorCode PetscOptionsAMSDestroy(void)
333d5649816SBarry Smith {
334d5649816SBarry Smith   PetscErrorCode ierr;
335d5649816SBarry Smith   AMS_Comm       acomm = -1;
336d5649816SBarry Smith   AMS_Memory     amem  = -1;
337d5649816SBarry Smith   char           options[16];
338d5649816SBarry Smith   const char     *string = "Exit";
339d5649816SBarry Smith 
340d5649816SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
341d5649816SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
342d5649816SBarry Smith   sprintf(options,"Options_%d",count++);
343d5649816SBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
344d5649816SBarry Smith   ierr = AMS_Memory_add_field(amem,"Exit",&string,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"Exit");
345d5649816SBarry Smith 
346d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
347d5649816SBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
348d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
349d5649816SBarry Smith   /* wait until accessor has unlocked the memory */
350d5649816SBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
351d5649816SBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
352d5649816SBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
353d5649816SBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
354d5649816SBarry Smith   PetscFunctionReturn(0);
355d5649816SBarry Smith }
356d5649816SBarry Smith 
357d5649816SBarry Smith #undef __FUNCT__
358d5649816SBarry Smith #define __FUNCT__ "PetscOptionsAMSInput"
359b3506946SBarry Smith /*
360d5649816SBarry Smith     PetscOptionsAMSInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the AMS
361b3506946SBarry Smith 
362b3506946SBarry Smith     Bugs:
363b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
364b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
365b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
366b3506946SBarry Smith 
367b3506946SBarry Smith 
368b3506946SBarry Smith */
369d5649816SBarry Smith PetscErrorCode PetscOptionsAMSInput()
370b3506946SBarry Smith {
371b3506946SBarry Smith   PetscErrorCode ierr;
372b3506946SBarry Smith   PetscOptions   next     = PetscOptionsObject.next;
373d5649816SBarry Smith   static int     mancount = 0;
374b3506946SBarry Smith   char           options[16];
375b3506946SBarry Smith   AMS_Comm       acomm         = -1;
376b3506946SBarry Smith   AMS_Memory     amem          = -1;
377ace3abfcSBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
3789f32e415SBarry Smith   char           manname[16];
379b3506946SBarry Smith 
380b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
381b3506946SBarry Smith   ierr = PetscViewerAMSGetAMSComm(PETSC_VIEWER_AMS_(PETSC_COMM_WORLD),&acomm);CHKERRQ(ierr);
382b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
3831ae3d29cSBarry Smith   ierr = AMS_Memory_create(acomm,options,&amem);CHKERRAMS(ierr);
3841ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
385a297a907SKarl Rupp 
3861bc75a8dSBarry Smith   PetscOptionsObject.pprefix = PetscOptionsObject.prefix; /* AMS will change this, so cannot pass prefix directly */
3871bc75a8dSBarry Smith 
3881ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,PetscOptionsObject.title,&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,PetscOptionsObject.title);
3890c74a584SJed Brown   /* ierr = AMS_Memory_add_field(amem,"mansec",&PetscOptionsObject.pprefix,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMS(ierr); */
3901ae3d29cSBarry Smith   ierr = AMS_Memory_add_field(amem,"ChangedMethod",&changedmethod,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,"ChangedMethod");
391b3506946SBarry Smith 
392b3506946SBarry Smith   while (next) {
3931ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,next->option,&next->set,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->option);
3941bc75a8dSBarry Smith     ierr =  PetscMalloc(sizeof(char*),&next->pman);CHKERRQ(ierr);
395a297a907SKarl Rupp 
3961bc75a8dSBarry Smith     *(char**)next->pman = next->man;
3979f32e415SBarry Smith     sprintf(manname,"man_%d",mancount++);
3981ae3d29cSBarry Smith     ierr = AMS_Memory_add_field(amem,manname,next->pman,1,AMS_STRING,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,manname);
3999f32e415SBarry Smith 
400b3506946SBarry Smith     switch (next->type) {
401b3506946SBarry Smith     case OPTION_HEAD:
402b3506946SBarry Smith       break;
403b3506946SBarry Smith     case OPTION_INT_ARRAY:
4041ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
405b3506946SBarry Smith       break;
406b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4071ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
408b3506946SBarry Smith       break;
409b3506946SBarry Smith     case OPTION_INT:
4101ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_INT,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
411b3506946SBarry Smith       break;
412b3506946SBarry Smith     case OPTION_REAL:
4131ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_DOUBLE,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
414b3506946SBarry Smith       break;
415b3506946SBarry Smith     case OPTION_LOGICAL:
4161ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4171ae3d29cSBarry Smith       break;
4181ae3d29cSBarry Smith     case OPTION_LOGICAL_ARRAY:
4191ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_BOOLEAN,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
42071f08665SBarry Smith       break;
421b3506946SBarry Smith     case OPTION_STRING:
4221ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4231ae3d29cSBarry Smith       break;
4241ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4251ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->data,next->arraylength,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
426b3506946SBarry Smith       break;
427b3506946SBarry Smith     case OPTION_LIST:
42871f08665SBarry Smith       {PetscInt ntext;
42971f08665SBarry Smith       char      ldefault[128];
43071f08665SBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
43171f08665SBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4321ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
433140e18c1SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
434d5649816SBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext-1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
4351ae3d29cSBarry Smith       break;}
4361ae3d29cSBarry Smith     case OPTION_ELIST:
437d5649816SBarry Smith       {PetscInt ntext = next->nlist;
4381ae3d29cSBarry Smith       char      ldefault[128];
4391ae3d29cSBarry Smith       ierr = PetscStrcpy(ldefault,"DEFAULT:");CHKERRQ(ierr);
4401ae3d29cSBarry Smith       ierr = PetscStrcat(ldefault,next->text);CHKERRQ(ierr);
4411ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,ldefault,next->data,1,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,ldefault);
4421ae3d29cSBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char**),&next->edata);CHKERRQ(ierr);
4431ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
4441ae3d29cSBarry Smith       ierr = AMS_Memory_add_field(amem,next->text,next->edata,ntext,AMS_STRING,AMS_WRITE,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRAMSFieldName(ierr,next->text);
44571f08665SBarry Smith       break;}
446b3506946SBarry Smith     default:
447b3506946SBarry Smith       break;
448b3506946SBarry Smith     }
449b3506946SBarry Smith     next = next->next;
450b3506946SBarry Smith   }
451b3506946SBarry Smith 
4521ae3d29cSBarry Smith   ierr = AMS_Memory_publish(amem);CHKERRAMS(ierr);
4531ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
454b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4551ae3d29cSBarry Smith   ierr = AMS_Memory_lock(amem,0);CHKERRAMS(ierr);
4561ae3d29cSBarry Smith   ierr = AMS_Memory_take_access(amem);CHKERRAMS(ierr);
457b3506946SBarry Smith 
458b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
459b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
460b3506946SBarry Smith 
4611ae3d29cSBarry Smith   ierr = AMS_Memory_grant_access(amem);CHKERRAMS(ierr);
4621ae3d29cSBarry Smith   ierr = AMS_Memory_destroy(amem);CHKERRAMS(ierr);
463b3506946SBarry Smith   PetscFunctionReturn(0);
464b3506946SBarry Smith }
465b3506946SBarry Smith #endif
466b3506946SBarry Smith 
4676356e834SBarry Smith #undef __FUNCT__
46853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
46953acd3b1SBarry Smith PetscErrorCode PetscOptionsEnd_Private(void)
47053acd3b1SBarry Smith {
47153acd3b1SBarry Smith   PetscErrorCode ierr;
4726356e834SBarry Smith   PetscOptions   last;
4736356e834SBarry Smith   char           option[256],value[1024],tmp[32];
474330cf3c9SBarry Smith   size_t         j;
47553acd3b1SBarry Smith 
47653acd3b1SBarry Smith   PetscFunctionBegin;
4771bc75a8dSBarry Smith   CHKMEMQ;
478aee2cecaSBarry Smith   if (PetscOptionsObject.next) {
479b3506946SBarry Smith     if (!PetscOptionsPublishCount) {
480b3506946SBarry Smith #if defined(PETSC_HAVE_AMS)
481d5649816SBarry Smith       ierr = PetscOptionsAMSInput();CHKERRQ(ierr);
482b3506946SBarry Smith #else
48371f08665SBarry Smith       ierr = PetscOptionsGetFromTextInput();CHKERRQ(ierr);
484b3506946SBarry Smith #endif
485aee2cecaSBarry Smith     }
486aee2cecaSBarry Smith   }
4876356e834SBarry Smith 
488c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.title);CHKERRQ(ierr);
489c31cb41cSBarry Smith   ierr = PetscFree(PetscOptionsObject.prefix);CHKERRQ(ierr);
4906356e834SBarry Smith 
4916356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
4926356e834SBarry Smith   if (PetscOptionsObject.changedmethod) PetscOptionsPublishCount = -2;
49361b37b28SSatish Balay   /* reset alreadyprinted flag */
49461b37b28SSatish Balay   PetscOptionsObject.alreadyprinted = PETSC_FALSE;
4953194b578SJed Brown   if (PetscOptionsObject.object) PetscOptionsObject.object->optionsprinted = PETSC_TRUE;
4960298fd71SBarry Smith   PetscOptionsObject.object = NULL;
4976356e834SBarry Smith 
4986356e834SBarry Smith   while (PetscOptionsObject.next) {
4996356e834SBarry Smith     if (PetscOptionsObject.next->set) {
5006356e834SBarry Smith       if (PetscOptionsObject.prefix) {
5016356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
5026356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.prefix);CHKERRQ(ierr);
5036356e834SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject.next->option+1);CHKERRQ(ierr);
5046356e834SBarry Smith       } else {
5056356e834SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject.next->option);CHKERRQ(ierr);
5066356e834SBarry Smith       }
5076356e834SBarry Smith 
5086356e834SBarry Smith       switch (PetscOptionsObject.next->type) {
5096356e834SBarry Smith       case OPTION_HEAD:
5106356e834SBarry Smith         break;
511e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
512e26ddf31SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[0]);
513e26ddf31SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
514e26ddf31SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject.next->data)[j]);
515e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
516e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
517e26ddf31SBarry Smith         }
518e26ddf31SBarry Smith         break;
5196356e834SBarry Smith       case OPTION_INT:
5207a72a596SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject.next->data);
5216356e834SBarry Smith         break;
5226356e834SBarry Smith       case OPTION_REAL:
5237a72a596SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject.next->data);
5246356e834SBarry Smith         break;
5256356e834SBarry Smith       case OPTION_REAL_ARRAY:
5267a72a596SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[0]);
5276356e834SBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5287a72a596SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject.next->data)[j]);
5296356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5306356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5316356e834SBarry Smith         }
5326356e834SBarry Smith         break;
5336356e834SBarry Smith       case OPTION_LOGICAL:
53471f08665SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject.next->data);
5356356e834SBarry Smith         break;
5361ae3d29cSBarry Smith       case OPTION_LOGICAL_ARRAY:
537ace3abfcSBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[0]);
5381ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
539ace3abfcSBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject.next->data)[j]);
5401ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5411ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5421ae3d29cSBarry Smith         }
5431ae3d29cSBarry Smith         break;
5446356e834SBarry Smith       case OPTION_LIST:
5451ae3d29cSBarry Smith       case OPTION_ELIST:
54671f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5476356e834SBarry Smith         break;
5481ae3d29cSBarry Smith       case OPTION_STRING:
54971f08665SBarry Smith         ierr = PetscStrcpy(value,*(char**)PetscOptionsObject.next->data);CHKERRQ(ierr);
5501ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
5511ae3d29cSBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject.next->data)[0]);
5521ae3d29cSBarry Smith         for (j=1; j<PetscOptionsObject.next->arraylength; j++) {
5531ae3d29cSBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject.next->data)[j]);
5541ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5551ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5561ae3d29cSBarry Smith         }
5576356e834SBarry Smith         break;
5586356e834SBarry Smith       }
5596356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5606356e834SBarry Smith     }
561503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->text);CHKERRQ(ierr);
562503cfb0cSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->option);CHKERRQ(ierr);
5636356e834SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->man);CHKERRQ(ierr);
56405b42c5fSBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->data);CHKERRQ(ierr);
56571f08665SBarry Smith     ierr   = PetscFree(PetscOptionsObject.next->edata);CHKERRQ(ierr);
566a297a907SKarl Rupp 
5676356e834SBarry Smith     last                    = PetscOptionsObject.next;
5686356e834SBarry Smith     PetscOptionsObject.next = PetscOptionsObject.next->next;
5696356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
5701bc75a8dSBarry Smith     CHKMEMQ;
5716356e834SBarry Smith   }
5721bc75a8dSBarry Smith   CHKMEMQ;
5736356e834SBarry Smith   PetscOptionsObject.next = 0;
57453acd3b1SBarry Smith   PetscFunctionReturn(0);
57553acd3b1SBarry Smith }
57653acd3b1SBarry Smith 
57753acd3b1SBarry Smith #undef __FUNCT__
57853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnum"
57953acd3b1SBarry Smith /*@C
58053acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
58153acd3b1SBarry Smith 
5823f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
58353acd3b1SBarry Smith 
58453acd3b1SBarry Smith    Input Parameters:
58553acd3b1SBarry Smith +  opt - option name
58653acd3b1SBarry Smith .  text - short string that describes the option
58753acd3b1SBarry Smith .  man - manual page with additional information on option
58853acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
58953acd3b1SBarry Smith -  defaultv - the default (current) value
59053acd3b1SBarry Smith 
59153acd3b1SBarry Smith    Output Parameter:
59253acd3b1SBarry Smith +  value - the  value to return
59353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
59453acd3b1SBarry Smith 
59553acd3b1SBarry Smith    Level: beginner
59653acd3b1SBarry Smith 
59753acd3b1SBarry Smith    Concepts: options database
59853acd3b1SBarry Smith 
59953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
60053acd3b1SBarry Smith 
60153acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
60253acd3b1SBarry Smith 
60353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
604acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
605acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
60653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
60753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
608acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
60953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
61053acd3b1SBarry Smith @*/
6117087cfbeSBarry Smith PetscErrorCode  PetscOptionsEnum(const char opt[],const char text[],const char man[],const char *const *list,PetscEnum defaultv,PetscEnum *value,PetscBool  *set)
61253acd3b1SBarry Smith {
61353acd3b1SBarry Smith   PetscErrorCode ierr;
61453acd3b1SBarry Smith   PetscInt       ntext = 0;
615aa5bb8c0SSatish Balay   PetscInt       tval;
616ace3abfcSBarry Smith   PetscBool      tflg;
61753acd3b1SBarry Smith 
61853acd3b1SBarry Smith   PetscFunctionBegin;
61953acd3b1SBarry Smith   while (list[ntext++]) {
620e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
62153acd3b1SBarry Smith   }
622e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
62353acd3b1SBarry Smith   ntext -= 3;
624aa5bb8c0SSatish Balay   ierr   = PetscOptionsEList(opt,text,man,list,ntext,list[defaultv],&tval,&tflg);CHKERRQ(ierr);
625aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
626aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
627aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
62853acd3b1SBarry Smith   PetscFunctionReturn(0);
62953acd3b1SBarry Smith }
63053acd3b1SBarry Smith 
63153acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
63253acd3b1SBarry Smith #undef __FUNCT__
63353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsInt"
63453acd3b1SBarry Smith /*@C
63553acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
63653acd3b1SBarry Smith 
6373f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
63853acd3b1SBarry Smith 
63953acd3b1SBarry Smith    Input Parameters:
64053acd3b1SBarry Smith +  opt - option name
64153acd3b1SBarry Smith .  text - short string that describes the option
64253acd3b1SBarry Smith .  man - manual page with additional information on option
64353acd3b1SBarry Smith -  defaultv - the default (current) value
64453acd3b1SBarry Smith 
64553acd3b1SBarry Smith    Output Parameter:
64653acd3b1SBarry Smith +  value - the integer value to return
64753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
64853acd3b1SBarry Smith 
64953acd3b1SBarry Smith    Level: beginner
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith    Concepts: options database^has int
65253acd3b1SBarry Smith 
65353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
65453acd3b1SBarry Smith 
65553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
656acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
657acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
65853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
65953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
660acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
66153acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
66253acd3b1SBarry Smith @*/
6637087cfbeSBarry Smith PetscErrorCode  PetscOptionsInt(const char opt[],const char text[],const char man[],PetscInt defaultv,PetscInt *value,PetscBool  *set)
66453acd3b1SBarry Smith {
66553acd3b1SBarry Smith   PetscErrorCode ierr;
6666356e834SBarry Smith   PetscOptions   amsopt;
66753acd3b1SBarry Smith 
66853acd3b1SBarry Smith   PetscFunctionBegin;
669b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
6706356e834SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
6716356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
672a297a907SKarl Rupp 
6736356e834SBarry Smith     *(PetscInt*)amsopt->data = defaultv;
674af6d86caSBarry Smith   }
67553acd3b1SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject.prefix,opt,value,set);CHKERRQ(ierr);
67661b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
6772aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,defaultv,text,ManSection(man));CHKERRQ(ierr);
67853acd3b1SBarry Smith   }
67953acd3b1SBarry Smith   PetscFunctionReturn(0);
68053acd3b1SBarry Smith }
68153acd3b1SBarry Smith 
68253acd3b1SBarry Smith #undef __FUNCT__
68353acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsString"
68453acd3b1SBarry Smith /*@C
68553acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
68653acd3b1SBarry Smith 
6873f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
68853acd3b1SBarry Smith 
68953acd3b1SBarry Smith    Input Parameters:
69053acd3b1SBarry Smith +  opt - option name
69153acd3b1SBarry Smith .  text - short string that describes the option
69253acd3b1SBarry Smith .  man - manual page with additional information on option
693bcbf2dc5SJed Brown .  defaultv - the default (current) value
694bcbf2dc5SJed Brown -  len - length of the result string including null terminator
69553acd3b1SBarry Smith 
69653acd3b1SBarry Smith    Output Parameter:
69753acd3b1SBarry Smith +  value - the value to return
69853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
69953acd3b1SBarry Smith 
70053acd3b1SBarry Smith    Level: beginner
70153acd3b1SBarry Smith 
70253acd3b1SBarry Smith    Concepts: options database^has int
70353acd3b1SBarry Smith 
70453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
70553acd3b1SBarry Smith 
7067fccdfe4SBarry 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).
7077fccdfe4SBarry Smith 
70853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
709acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
710acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
71153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
71253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
713acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
71453acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
71553acd3b1SBarry Smith @*/
7167087cfbeSBarry Smith PetscErrorCode  PetscOptionsString(const char opt[],const char text[],const char man[],const char defaultv[],char value[],size_t len,PetscBool  *set)
71753acd3b1SBarry Smith {
71853acd3b1SBarry Smith   PetscErrorCode ierr;
719aee2cecaSBarry Smith   PetscOptions   amsopt;
72053acd3b1SBarry Smith 
72153acd3b1SBarry Smith   PetscFunctionBegin;
722b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
723aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
72471f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
725a297a907SKarl Rupp 
72671f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
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(),
76453acd3b1SBarry Smith           PetscOptionsList(), 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(),
81453acd3b1SBarry Smith           PetscOptionsList(), 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(),
85753acd3b1SBarry Smith           PetscOptionsList(), 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) {
8661ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&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__
87953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsList"
88053acd3b1SBarry Smith /*@C
88153acd3b1SBarry Smith      PetscOptionsList - 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 
90653acd3b1SBarry Smith    Concepts: options database^list
90753acd3b1SBarry Smith 
90853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
909acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
91053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
91153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
912acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
913171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsEnum()
91453acd3b1SBarry Smith @*/
915140e18c1SBarry Smith PetscErrorCode  PetscOptionsList(const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char defaultv[],char value[],size_t len,PetscBool  *set)
91653acd3b1SBarry Smith {
91753acd3b1SBarry Smith   PetscErrorCode ierr;
9183cc1e11dSBarry Smith   PetscOptions   amsopt;
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith   PetscFunctionBegin;
921b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
9223cc1e11dSBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_LIST,&amsopt);CHKERRQ(ierr);
92371f08665SBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
924a297a907SKarl Rupp 
92571f08665SBarry Smith     *(const char**)amsopt->data = defaultv;
926a297a907SKarl Rupp 
9273cc1e11dSBarry Smith     amsopt->flist = list;
9283cc1e11dSBarry Smith   }
92953acd3b1SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject.prefix,opt,value,len,set);CHKERRQ(ierr);
93061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
931140e18c1SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject.comm,stdout,PetscOptionsObject.prefix,opt,ltext,man,list,defaultv);CHKERRQ(ierr);CHKERRQ(ierr);
93253acd3b1SBarry Smith   }
93353acd3b1SBarry Smith   PetscFunctionReturn(0);
93453acd3b1SBarry Smith }
93553acd3b1SBarry Smith 
93653acd3b1SBarry Smith #undef __FUNCT__
93753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEList"
93853acd3b1SBarry Smith /*@C
93953acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
94053acd3b1SBarry Smith 
9413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
94253acd3b1SBarry Smith 
94353acd3b1SBarry Smith    Input Parameters:
94453acd3b1SBarry Smith +  opt - option name
94553acd3b1SBarry Smith .  ltext - short string that describes the option
94653acd3b1SBarry Smith .  man - manual page with additional information on option
94753acd3b1SBarry Smith .  list - the possible choices
94853acd3b1SBarry Smith .  ntext - number of choices
94953acd3b1SBarry Smith -  defaultv - the default (current) value
95053acd3b1SBarry Smith 
95153acd3b1SBarry Smith    Output Parameter:
95253acd3b1SBarry Smith +  value - the index of the value to return
95353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
95453acd3b1SBarry Smith 
95553acd3b1SBarry Smith    Level: intermediate
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
95853acd3b1SBarry Smith 
959140e18c1SBarry Smith    See PetscOptionsList() for when the choices are given in a PetscFunctionList()
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith    Concepts: options database^list
96253acd3b1SBarry Smith 
96353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
964acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
96553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
967acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
968171400e9SBarry Smith           PetscOptionsList(), PetscOptionsEnum()
96953acd3b1SBarry Smith @*/
9707087cfbeSBarry 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)
97153acd3b1SBarry Smith {
97253acd3b1SBarry Smith   PetscErrorCode ierr;
97353acd3b1SBarry Smith   PetscInt       i;
9741ae3d29cSBarry Smith   PetscOptions   amsopt;
97553acd3b1SBarry Smith 
97653acd3b1SBarry Smith   PetscFunctionBegin;
9771ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
978d5649816SBarry Smith     ierr = PetscOptionsCreate_Private(opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
9791ae3d29cSBarry Smith     ierr = PetscMalloc(sizeof(char*),&amsopt->data);CHKERRQ(ierr);
980a297a907SKarl Rupp 
9811ae3d29cSBarry Smith     *(const char**)amsopt->data = defaultv;
982a297a907SKarl Rupp 
9831ae3d29cSBarry Smith     amsopt->list  = list;
9841ae3d29cSBarry Smith     amsopt->nlist = ntext;
9851ae3d29cSBarry Smith   }
98653acd3b1SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject.prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
98761b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
98853acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,defaultv);CHKERRQ(ierr);
98953acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
99053acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," %s",list[i]);CHKERRQ(ierr);
99153acd3b1SBarry Smith     }
9922aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
99353acd3b1SBarry Smith   }
99453acd3b1SBarry Smith   PetscFunctionReturn(0);
99553acd3b1SBarry Smith }
99653acd3b1SBarry Smith 
99753acd3b1SBarry Smith #undef __FUNCT__
998acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupBegin"
99953acd3b1SBarry Smith /*@C
1000acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1001d5649816SBarry Smith        which at most a single value can be true.
100253acd3b1SBarry Smith 
10033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Input Parameters:
100653acd3b1SBarry Smith +  opt - option name
100753acd3b1SBarry Smith .  text - short string that describes the option
100853acd3b1SBarry Smith -  man - manual page with additional information on option
100953acd3b1SBarry Smith 
101053acd3b1SBarry Smith    Output Parameter:
101153acd3b1SBarry Smith .  flg - whether that option was set or not
101253acd3b1SBarry Smith 
101353acd3b1SBarry Smith    Level: intermediate
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
101653acd3b1SBarry Smith 
1017acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
101853acd3b1SBarry Smith 
101953acd3b1SBarry Smith     Concepts: options database^logical group
102053acd3b1SBarry Smith 
102153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1022acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
102353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
102453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1025acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
102653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
102753acd3b1SBarry Smith @*/
10287087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin(const char opt[],const char text[],const char man[],PetscBool  *flg)
102953acd3b1SBarry Smith {
103053acd3b1SBarry Smith   PetscErrorCode ierr;
10311ae3d29cSBarry Smith   PetscOptions   amsopt;
103253acd3b1SBarry Smith 
103353acd3b1SBarry Smith   PetscFunctionBegin;
10341ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10351ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1036ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1037a297a907SKarl Rupp 
1038ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10391ae3d29cSBarry Smith   }
104068b16fdaSBarry Smith   *flg = PETSC_FALSE;
10410298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
104261b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
104353acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
10442aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
104553acd3b1SBarry Smith   }
104653acd3b1SBarry Smith   PetscFunctionReturn(0);
104753acd3b1SBarry Smith }
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith #undef __FUNCT__
1050acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroup"
105153acd3b1SBarry Smith /*@C
1052acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1053d5649816SBarry Smith        which at most a single value can be true.
105453acd3b1SBarry Smith 
10553f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105653acd3b1SBarry Smith 
105753acd3b1SBarry Smith    Input Parameters:
105853acd3b1SBarry Smith +  opt - option name
105953acd3b1SBarry Smith .  text - short string that describes the option
106053acd3b1SBarry Smith -  man - manual page with additional information on option
106153acd3b1SBarry Smith 
106253acd3b1SBarry Smith    Output Parameter:
106353acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
106453acd3b1SBarry Smith 
106553acd3b1SBarry Smith    Level: intermediate
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106853acd3b1SBarry Smith 
1069acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
107053acd3b1SBarry Smith 
107153acd3b1SBarry Smith     Concepts: options database^logical group
107253acd3b1SBarry Smith 
107353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1074acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1077acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
107853acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
107953acd3b1SBarry Smith @*/
10807087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroup(const char opt[],const char text[],const char man[],PetscBool  *flg)
108153acd3b1SBarry Smith {
108253acd3b1SBarry Smith   PetscErrorCode ierr;
10831ae3d29cSBarry Smith   PetscOptions   amsopt;
108453acd3b1SBarry Smith 
108553acd3b1SBarry Smith   PetscFunctionBegin;
10861ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
10871ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1088ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1089a297a907SKarl Rupp 
1090ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10911ae3d29cSBarry Smith   }
109217326d04SJed Brown   *flg = PETSC_FALSE;
10930298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
109461b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
10952aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109653acd3b1SBarry Smith   }
109753acd3b1SBarry Smith   PetscFunctionReturn(0);
109853acd3b1SBarry Smith }
109953acd3b1SBarry Smith 
110053acd3b1SBarry Smith #undef __FUNCT__
1101acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolGroupEnd"
110253acd3b1SBarry Smith /*@C
1103acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1104d5649816SBarry Smith        which at most a single value can be true.
110553acd3b1SBarry Smith 
11063f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110753acd3b1SBarry Smith 
110853acd3b1SBarry Smith    Input Parameters:
110953acd3b1SBarry Smith +  opt - option name
111053acd3b1SBarry Smith .  text - short string that describes the option
111153acd3b1SBarry Smith -  man - manual page with additional information on option
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Output Parameter:
111453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111553acd3b1SBarry Smith 
111653acd3b1SBarry Smith    Level: intermediate
111753acd3b1SBarry Smith 
111853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111953acd3b1SBarry Smith 
1120acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
112153acd3b1SBarry Smith 
112253acd3b1SBarry Smith     Concepts: options database^logical group
112353acd3b1SBarry Smith 
112453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1125acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1128acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
112953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
113053acd3b1SBarry Smith @*/
11317087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd(const char opt[],const char text[],const char man[],PetscBool  *flg)
113253acd3b1SBarry Smith {
113353acd3b1SBarry Smith   PetscErrorCode ierr;
11341ae3d29cSBarry Smith   PetscOptions   amsopt;
113553acd3b1SBarry Smith 
113653acd3b1SBarry Smith   PetscFunctionBegin;
11371ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
11381ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1139ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1140a297a907SKarl Rupp 
1141ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11421ae3d29cSBarry Smith   }
114317326d04SJed Brown   *flg = PETSC_FALSE;
11440298fd71SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,NULL);CHKERRQ(ierr);
114561b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
11462aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"    -%s%s: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114753acd3b1SBarry Smith   }
114853acd3b1SBarry Smith   PetscFunctionReturn(0);
114953acd3b1SBarry Smith }
115053acd3b1SBarry Smith 
115153acd3b1SBarry Smith #undef __FUNCT__
1152acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBool"
115353acd3b1SBarry Smith /*@C
1154acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
115553acd3b1SBarry Smith 
11563f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115753acd3b1SBarry Smith 
115853acd3b1SBarry Smith    Input Parameters:
115953acd3b1SBarry Smith +  opt - option name
116053acd3b1SBarry Smith .  text - short string that describes the option
116153acd3b1SBarry Smith -  man - manual page with additional information on option
116253acd3b1SBarry Smith 
116353acd3b1SBarry Smith    Output Parameter:
116453acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
116553acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith    Level: beginner
116853acd3b1SBarry Smith 
116953acd3b1SBarry Smith    Concepts: options database^logical
117053acd3b1SBarry Smith 
117153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
117253acd3b1SBarry Smith 
117353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1174acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1175acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
117653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1178acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
117953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
118053acd3b1SBarry Smith @*/
11817087cfbeSBarry Smith PetscErrorCode  PetscOptionsBool(const char opt[],const char text[],const char man[],PetscBool deflt,PetscBool  *flg,PetscBool  *set)
118253acd3b1SBarry Smith {
118353acd3b1SBarry Smith   PetscErrorCode ierr;
1184ace3abfcSBarry Smith   PetscBool      iset;
1185aee2cecaSBarry Smith   PetscOptions   amsopt;
118653acd3b1SBarry Smith 
118753acd3b1SBarry Smith   PetscFunctionBegin;
1188b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1189aee2cecaSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL,&amsopt);CHKERRQ(ierr);
1190ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1191a297a907SKarl Rupp 
1192ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = deflt;
1193af6d86caSBarry Smith   }
1194acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PetscOptionsObject.prefix,opt,flg,&iset);CHKERRQ(ierr);
119553acd3b1SBarry Smith   if (!iset) {
119653acd3b1SBarry Smith     if (flg) *flg = deflt;
119753acd3b1SBarry Smith   }
119853acd3b1SBarry Smith   if (set) *set = iset;
119961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1200ace3abfcSBarry Smith     const char *v = PetscBools[deflt];
12012aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
120253acd3b1SBarry Smith   }
120353acd3b1SBarry Smith   PetscFunctionReturn(0);
120453acd3b1SBarry Smith }
120553acd3b1SBarry Smith 
120653acd3b1SBarry Smith #undef __FUNCT__
120753acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsRealArray"
120853acd3b1SBarry Smith /*@C
120953acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
121053acd3b1SBarry Smith    option in the database. The values must be separated with commas with
121153acd3b1SBarry Smith    no intervening spaces.
121253acd3b1SBarry Smith 
12133f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
121453acd3b1SBarry Smith 
121553acd3b1SBarry Smith    Input Parameters:
121653acd3b1SBarry Smith +  opt - the option one is seeking
121753acd3b1SBarry Smith .  text - short string describing option
121853acd3b1SBarry Smith .  man - manual page for option
121953acd3b1SBarry Smith -  nmax - maximum number of values
122053acd3b1SBarry Smith 
122153acd3b1SBarry Smith    Output Parameter:
122253acd3b1SBarry Smith +  value - location to copy values
122353acd3b1SBarry Smith .  nmax - actual number of values found
122453acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
122553acd3b1SBarry Smith 
122653acd3b1SBarry Smith    Level: beginner
122753acd3b1SBarry Smith 
122853acd3b1SBarry Smith    Notes:
122953acd3b1SBarry Smith    The user should pass in an array of doubles
123053acd3b1SBarry Smith 
123153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
123253acd3b1SBarry Smith 
123353acd3b1SBarry Smith    Concepts: options database^array of strings
123453acd3b1SBarry Smith 
123553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1236acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
123753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
123853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1239acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
124053acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
124153acd3b1SBarry Smith @*/
12427087cfbeSBarry Smith PetscErrorCode  PetscOptionsRealArray(const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
124353acd3b1SBarry Smith {
124453acd3b1SBarry Smith   PetscErrorCode ierr;
124553acd3b1SBarry Smith   PetscInt       i;
1246e26ddf31SBarry Smith   PetscOptions   amsopt;
124753acd3b1SBarry Smith 
124853acd3b1SBarry Smith   PetscFunctionBegin;
1249b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1250e26ddf31SBarry Smith     PetscReal *vals;
1251e26ddf31SBarry Smith 
1252e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1253e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1254e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1255e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1256e26ddf31SBarry Smith     amsopt->arraylength = *n;
1257e26ddf31SBarry Smith   }
125853acd3b1SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
125961b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1260a83599f4SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%G",PetscOptionsObject.prefix?PetscOptionsObject.prefix:"",opt+1,value[0]);CHKERRQ(ierr);
126153acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1262a83599f4SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%G",value[i]);CHKERRQ(ierr);
126353acd3b1SBarry Smith     }
12642aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
126553acd3b1SBarry Smith   }
126653acd3b1SBarry Smith   PetscFunctionReturn(0);
126753acd3b1SBarry Smith }
126853acd3b1SBarry Smith 
126953acd3b1SBarry Smith 
127053acd3b1SBarry Smith #undef __FUNCT__
127153acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsIntArray"
127253acd3b1SBarry Smith /*@C
127353acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1274b32a342fSShri Abhyankar    option in the database.
127553acd3b1SBarry Smith 
12763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
127753acd3b1SBarry Smith 
127853acd3b1SBarry Smith    Input Parameters:
127953acd3b1SBarry Smith +  opt - the option one is seeking
128053acd3b1SBarry Smith .  text - short string describing option
128153acd3b1SBarry Smith .  man - manual page for option
1282f8a50e2bSBarry Smith -  n - maximum number of values
128353acd3b1SBarry Smith 
128453acd3b1SBarry Smith    Output Parameter:
128553acd3b1SBarry Smith +  value - location to copy values
1286f8a50e2bSBarry Smith .  n - actual number of values found
128753acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
128853acd3b1SBarry Smith 
128953acd3b1SBarry Smith    Level: beginner
129053acd3b1SBarry Smith 
129153acd3b1SBarry Smith    Notes:
1292b32a342fSShri Abhyankar    The array can be passed as
1293b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
12940fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
12950fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
12960fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1297b32a342fSShri Abhyankar 
1298b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
129953acd3b1SBarry Smith 
130053acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
130153acd3b1SBarry Smith 
1302b32a342fSShri Abhyankar    Concepts: options database^array of ints
130353acd3b1SBarry Smith 
130453acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1305acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
130653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
130753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1308acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
130953acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList(), PetscOptionsRealArray()
131053acd3b1SBarry Smith @*/
13117087cfbeSBarry Smith PetscErrorCode  PetscOptionsIntArray(const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
131253acd3b1SBarry Smith {
131353acd3b1SBarry Smith   PetscErrorCode ierr;
131453acd3b1SBarry Smith   PetscInt       i;
1315e26ddf31SBarry Smith   PetscOptions   amsopt;
131653acd3b1SBarry Smith 
131753acd3b1SBarry Smith   PetscFunctionBegin;
1318b3506946SBarry Smith   if (!PetscOptionsPublishCount) {
1319e26ddf31SBarry Smith     PetscInt *vals;
1320e26ddf31SBarry Smith 
1321e26ddf31SBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1322e26ddf31SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1323e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1324e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1325e26ddf31SBarry Smith     amsopt->arraylength = *n;
1326e26ddf31SBarry Smith   }
132753acd3b1SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
132861b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
132953acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
133053acd3b1SBarry Smith     for (i=1; i<*n; i++) {
133153acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
133253acd3b1SBarry Smith     }
13332aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
133453acd3b1SBarry Smith   }
133553acd3b1SBarry Smith   PetscFunctionReturn(0);
133653acd3b1SBarry Smith }
133753acd3b1SBarry Smith 
133853acd3b1SBarry Smith #undef __FUNCT__
133953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsStringArray"
134053acd3b1SBarry Smith /*@C
134153acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
134253acd3b1SBarry Smith    option in the database. The values must be separated with commas with
134353acd3b1SBarry Smith    no intervening spaces.
134453acd3b1SBarry Smith 
13453f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
134653acd3b1SBarry Smith 
134753acd3b1SBarry Smith    Input Parameters:
134853acd3b1SBarry Smith +  opt - the option one is seeking
134953acd3b1SBarry Smith .  text - short string describing option
135053acd3b1SBarry Smith .  man - manual page for option
135153acd3b1SBarry Smith -  nmax - maximum number of strings
135253acd3b1SBarry Smith 
135353acd3b1SBarry Smith    Output Parameter:
135453acd3b1SBarry Smith +  value - location to copy strings
135553acd3b1SBarry Smith .  nmax - actual number of strings found
135653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
135753acd3b1SBarry Smith 
135853acd3b1SBarry Smith    Level: beginner
135953acd3b1SBarry Smith 
136053acd3b1SBarry Smith    Notes:
136153acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
136253acd3b1SBarry Smith    strings returned by this function.
136353acd3b1SBarry Smith 
136453acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
136553acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
136653acd3b1SBarry Smith 
136753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
136853acd3b1SBarry Smith 
136953acd3b1SBarry Smith    Concepts: options database^array of strings
137053acd3b1SBarry Smith 
137153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1372acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
137353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
137453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1375acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
137653acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
137753acd3b1SBarry Smith @*/
13787087cfbeSBarry Smith PetscErrorCode  PetscOptionsStringArray(const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
137953acd3b1SBarry Smith {
138053acd3b1SBarry Smith   PetscErrorCode ierr;
13811ae3d29cSBarry Smith   PetscOptions   amsopt;
138253acd3b1SBarry Smith 
138353acd3b1SBarry Smith   PetscFunctionBegin;
13841ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
13851ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
13861ae3d29cSBarry Smith     ierr = PetscMalloc((*nmax)*sizeof(char*),&amsopt->data);CHKERRQ(ierr);
1387a297a907SKarl Rupp 
13881ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
13891ae3d29cSBarry Smith   }
139053acd3b1SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject.prefix,opt,value,nmax,set);CHKERRQ(ierr);
139161b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
13922aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
139353acd3b1SBarry Smith   }
139453acd3b1SBarry Smith   PetscFunctionReturn(0);
139553acd3b1SBarry Smith }
139653acd3b1SBarry Smith 
1397e2446a98SMatthew Knepley #undef __FUNCT__
1398acfcf0e5SJed Brown #define __FUNCT__ "PetscOptionsBoolArray"
1399e2446a98SMatthew Knepley /*@C
1400acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1401e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1402e2446a98SMatthew Knepley    no intervening spaces.
1403e2446a98SMatthew Knepley 
14043f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1405e2446a98SMatthew Knepley 
1406e2446a98SMatthew Knepley    Input Parameters:
1407e2446a98SMatthew Knepley +  opt - the option one is seeking
1408e2446a98SMatthew Knepley .  text - short string describing option
1409e2446a98SMatthew Knepley .  man - manual page for option
1410e2446a98SMatthew Knepley -  nmax - maximum number of values
1411e2446a98SMatthew Knepley 
1412e2446a98SMatthew Knepley    Output Parameter:
1413e2446a98SMatthew Knepley +  value - location to copy values
1414e2446a98SMatthew Knepley .  nmax - actual number of values found
1415e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1416e2446a98SMatthew Knepley 
1417e2446a98SMatthew Knepley    Level: beginner
1418e2446a98SMatthew Knepley 
1419e2446a98SMatthew Knepley    Notes:
1420e2446a98SMatthew Knepley    The user should pass in an array of doubles
1421e2446a98SMatthew Knepley 
1422e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1423e2446a98SMatthew Knepley 
1424e2446a98SMatthew Knepley    Concepts: options database^array of strings
1425e2446a98SMatthew Knepley 
1426e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1427acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1428e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1429e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1430acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1431e2446a98SMatthew Knepley           PetscOptionsList(), PetscOptionsEList()
1432e2446a98SMatthew Knepley @*/
14337087cfbeSBarry Smith PetscErrorCode  PetscOptionsBoolArray(const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1434e2446a98SMatthew Knepley {
1435e2446a98SMatthew Knepley   PetscErrorCode ierr;
1436e2446a98SMatthew Knepley   PetscInt       i;
14371ae3d29cSBarry Smith   PetscOptions   amsopt;
1438e2446a98SMatthew Knepley 
1439e2446a98SMatthew Knepley   PetscFunctionBegin;
14401ae3d29cSBarry Smith   if (!PetscOptionsPublishCount) {
1441ace3abfcSBarry Smith     PetscBool *vals;
14421ae3d29cSBarry Smith 
14431ae3d29cSBarry Smith     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_LOGICAL_ARRAY,&amsopt);CHKERRQ(ierr);
1444ace3abfcSBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1445ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14461ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14471ae3d29cSBarry Smith     amsopt->arraylength = *n;
14481ae3d29cSBarry Smith   }
1449acfcf0e5SJed Brown   ierr = PetscOptionsGetBoolArray(PetscOptionsObject.prefix,opt,value,n,set);CHKERRQ(ierr);
1450e2446a98SMatthew Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
1451e2446a98SMatthew Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%d",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1452e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1453e2446a98SMatthew Knepley       ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,",%d",value[i]);CHKERRQ(ierr);
1454e2446a98SMatthew Knepley     }
14552aa6d131SJed Brown     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1456e2446a98SMatthew Knepley   }
1457e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1458e2446a98SMatthew Knepley }
1459e2446a98SMatthew Knepley 
14608cc676e6SMatthew G Knepley #undef __FUNCT__
14618cc676e6SMatthew G Knepley #define __FUNCT__ "PetscOptionsViewer"
14628cc676e6SMatthew G Knepley /*@C
14638cc676e6SMatthew G Knepley    PetscOptionsInt - Gets a viewer appropriate for the type indicated by the user
14648cc676e6SMatthew G Knepley 
14658cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
14668cc676e6SMatthew G Knepley 
14678cc676e6SMatthew G Knepley    Input Parameters:
14688cc676e6SMatthew G Knepley +  opt - option name
14698cc676e6SMatthew G Knepley .  text - short string that describes the option
14708cc676e6SMatthew G Knepley -  man - manual page with additional information on option
14718cc676e6SMatthew G Knepley 
14728cc676e6SMatthew G Knepley    Output Parameter:
14738cc676e6SMatthew G Knepley +  viewer - the viewer
14748cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
14758cc676e6SMatthew G Knepley 
14768cc676e6SMatthew G Knepley    Level: beginner
14778cc676e6SMatthew G Knepley 
14788cc676e6SMatthew G Knepley    Concepts: options database^has int
14798cc676e6SMatthew G Knepley 
14808cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
14818cc676e6SMatthew G Knepley If no value is provided ascii:stdout is used
14828cc676e6SMatthew 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
14838cc676e6SMatthew G Knepley $                                     about the object to standard out
14848cc676e6SMatthew G Knepley $       binary[:filename]   defaults to binaryoutput
14858cc676e6SMatthew G Knepley $       draw
14868cc676e6SMatthew G Knepley $       socket[:port]    defaults to the standard output port
14878cc676e6SMatthew G Knepley 
1488cffb1e40SBarry Smith    Use PetscRestoreViewerDestroy() after using the viewer, otherwise a memory leak will occur
14898cc676e6SMatthew G Knepley 
14908cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
14918cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
14928cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
14938cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
14948cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
14958cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
14968cc676e6SMatthew G Knepley           PetscOptionsList(), PetscOptionsEList()
14978cc676e6SMatthew G Knepley @*/
1498cffb1e40SBarry Smith PetscErrorCode  PetscOptionsViewer(const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
14998cc676e6SMatthew G Knepley {
15008cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15018cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15028cc676e6SMatthew G Knepley 
15038cc676e6SMatthew G Knepley   PetscFunctionBegin;
15048cc676e6SMatthew G Knepley   if (!PetscOptionsPublishCount) {
15058cc676e6SMatthew G Knepley     ierr = PetscOptionsCreate_Private(opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
15068cc676e6SMatthew G Knepley     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
1507a297a907SKarl Rupp 
15088cc676e6SMatthew G Knepley     *(const char**)amsopt->data = "";
15098cc676e6SMatthew G Knepley   }
1510cffb1e40SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject.comm,PetscOptionsObject.prefix,opt,viewer,format,set);CHKERRQ(ierr);
15118cc676e6SMatthew G Knepley   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
15128cc676e6SMatthew G Knepley     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject.prefix ? PetscOptionsObject.prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15138cc676e6SMatthew G Knepley   }
15148cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15158cc676e6SMatthew G Knepley }
15168cc676e6SMatthew G Knepley 
151753acd3b1SBarry Smith 
151853acd3b1SBarry Smith #undef __FUNCT__
151953acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
152053acd3b1SBarry Smith /*@C
1521b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
152253acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
152353acd3b1SBarry Smith 
15243f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
152553acd3b1SBarry Smith 
152653acd3b1SBarry Smith    Input Parameter:
152753acd3b1SBarry Smith .   head - the heading text
152853acd3b1SBarry Smith 
152953acd3b1SBarry Smith 
153053acd3b1SBarry Smith    Level: intermediate
153153acd3b1SBarry Smith 
153253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
153353acd3b1SBarry Smith 
1534b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
153553acd3b1SBarry Smith 
153653acd3b1SBarry Smith    Concepts: options database^subheading
153753acd3b1SBarry Smith 
153853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1539acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
154053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
154153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1542acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
154353acd3b1SBarry Smith           PetscOptionsList(), PetscOptionsEList()
154453acd3b1SBarry Smith @*/
15457087cfbeSBarry Smith PetscErrorCode  PetscOptionsHead(const char head[])
154653acd3b1SBarry Smith {
154753acd3b1SBarry Smith   PetscErrorCode ierr;
154853acd3b1SBarry Smith 
154953acd3b1SBarry Smith   PetscFunctionBegin;
155061b37b28SSatish Balay   if (PetscOptionsObject.printhelp && PetscOptionsPublishCount == 1 && !PetscOptionsObject.alreadyprinted) {
155153acd3b1SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject.comm,"  %s\n",head);CHKERRQ(ierr);
155253acd3b1SBarry Smith   }
155353acd3b1SBarry Smith   PetscFunctionReturn(0);
155453acd3b1SBarry Smith }
155553acd3b1SBarry Smith 
155653acd3b1SBarry Smith 
155753acd3b1SBarry Smith 
155853acd3b1SBarry Smith 
155953acd3b1SBarry Smith 
156053acd3b1SBarry Smith 
1561