xref: /petsc/src/sys/objects/aoptions.c (revision 050cccc3bd2925a9762563a38f4e695b2d452c28)
17d0a6c19SBarry Smith 
21a1499c8SBarry Smith 
353acd3b1SBarry Smith /*
43fc1eb6aSBarry Smith    Implements the higher-level options database querying methods. These are self-documenting and can attach at runtime to
53fc1eb6aSBarry Smith    GUI code to display the options and get values from the users.
653acd3b1SBarry Smith 
753acd3b1SBarry Smith */
853acd3b1SBarry Smith 
9afcb2eb5SJed Brown #include <petsc-private/petscimpl.h>        /*I  "petscsys.h"   I*/
10665c2dedSJed Brown #include <petscviewer.h>
1153acd3b1SBarry Smith 
122aa6d131SJed Brown #define ManSection(str) ((str) ? (str) : "None")
132aa6d131SJed Brown 
1453acd3b1SBarry Smith /*
1553acd3b1SBarry Smith     Keep a linked list of options that have been posted and we are waiting for
163fc1eb6aSBarry Smith    user selection. See the manual page for PetscOptionsBegin()
1753acd3b1SBarry Smith 
1853acd3b1SBarry Smith     Eventually we'll attach this beast to a MPI_Comm
1953acd3b1SBarry Smith */
20e55864a3SBarry Smith 
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 */
278c34d3f5SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptions *PetscOptionsObject,MPI_Comm comm,const char prefix[],const char title[],const char mansec[])
2853acd3b1SBarry Smith {
2953acd3b1SBarry Smith   PetscErrorCode ierr;
3053acd3b1SBarry Smith 
3153acd3b1SBarry Smith   PetscFunctionBegin;
32e55864a3SBarry Smith   PetscOptionsObject->next          = 0;
33e55864a3SBarry Smith   PetscOptionsObject->comm          = comm;
34e55864a3SBarry Smith   PetscOptionsObject->changedmethod = PETSC_FALSE;
35a297a907SKarl Rupp 
36e55864a3SBarry Smith   ierr = PetscStrallocpy(prefix,&PetscOptionsObject->prefix);CHKERRQ(ierr);
37e55864a3SBarry Smith   ierr = PetscStrallocpy(title,&PetscOptionsObject->title);CHKERRQ(ierr);
3853acd3b1SBarry Smith 
39e55864a3SBarry Smith   ierr = PetscOptionsHasName(NULL,"-help",&PetscOptionsObject->printhelp);CHKERRQ(ierr);
40e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1) {
41e55864a3SBarry Smith     if (!PetscOptionsObject->alreadyprinted) {
4253acd3b1SBarry Smith       ierr = (*PetscHelpPrintf)(comm,"%s -------------------------------------------------\n",title);CHKERRQ(ierr);
4353acd3b1SBarry Smith     }
4461b37b28SSatish Balay   }
4553acd3b1SBarry Smith   PetscFunctionReturn(0);
4653acd3b1SBarry Smith }
4753acd3b1SBarry Smith 
483194b578SJed Brown #undef __FUNCT__
493194b578SJed Brown #define __FUNCT__ "PetscObjectOptionsBegin_Private"
503194b578SJed Brown /*
513194b578SJed Brown     Handles setting up the data structure in a call to PetscObjectOptionsBegin()
523194b578SJed Brown */
538c34d3f5SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptions *PetscOptionsObject,PetscObject obj)
543194b578SJed Brown {
553194b578SJed Brown   PetscErrorCode ierr;
563194b578SJed Brown   char           title[256];
573194b578SJed Brown   PetscBool      flg;
583194b578SJed Brown 
593194b578SJed Brown   PetscFunctionBegin;
603194b578SJed Brown   PetscValidHeader(obj,1);
61e55864a3SBarry Smith   PetscOptionsObject->object         = obj;
62e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = obj->optionsprinted;
63a297a907SKarl Rupp 
643194b578SJed Brown   ierr = PetscStrcmp(obj->description,obj->class_name,&flg);CHKERRQ(ierr);
653194b578SJed Brown   if (flg) {
668caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s options",obj->class_name);CHKERRQ(ierr);
673194b578SJed Brown   } else {
688caf3d72SBarry Smith     ierr = PetscSNPrintf(title,sizeof(title),"%s (%s) options",obj->description,obj->class_name);CHKERRQ(ierr);
693194b578SJed Brown   }
70e55864a3SBarry Smith   ierr = PetscOptionsBegin_Private(PetscOptionsObject,obj->comm,obj->prefix,title,obj->mansec);CHKERRQ(ierr);
713194b578SJed Brown   PetscFunctionReturn(0);
723194b578SJed Brown }
733194b578SJed Brown 
7453acd3b1SBarry Smith /*
7553acd3b1SBarry Smith      Handles adding another option to the list of options within this particular PetscOptionsBegin() PetscOptionsEnd()
7653acd3b1SBarry Smith */
7753acd3b1SBarry Smith #undef __FUNCT__
7853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsCreate_Private"
798c34d3f5SBarry Smith static int PetscOptionsCreate_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOption *amsopt)
8053acd3b1SBarry Smith {
8153acd3b1SBarry Smith   int          ierr;
828c34d3f5SBarry Smith   PetscOption next;
833be6e4c3SJed Brown   PetscBool    valid;
8453acd3b1SBarry Smith 
8553acd3b1SBarry Smith   PetscFunctionBegin;
863be6e4c3SJed Brown   ierr = PetscOptionsValidKey(opt,&valid);CHKERRQ(ierr);
873be6e4c3SJed Brown   if (!valid) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"The option '%s' is not a valid key",opt);
883be6e4c3SJed Brown 
89b00a9115SJed Brown   ierr            = PetscNew(amsopt);CHKERRQ(ierr);
9053acd3b1SBarry Smith   (*amsopt)->next = 0;
9153acd3b1SBarry Smith   (*amsopt)->set  = PETSC_FALSE;
926356e834SBarry Smith   (*amsopt)->type = t;
9353acd3b1SBarry Smith   (*amsopt)->data = 0;
9461b37b28SSatish Balay 
9553acd3b1SBarry Smith   ierr = PetscStrallocpy(text,&(*amsopt)->text);CHKERRQ(ierr);
9653acd3b1SBarry Smith   ierr = PetscStrallocpy(opt,&(*amsopt)->option);CHKERRQ(ierr);
976356e834SBarry Smith   ierr = PetscStrallocpy(man,&(*amsopt)->man);CHKERRQ(ierr);
9853acd3b1SBarry Smith 
99e55864a3SBarry Smith   if (!PetscOptionsObject->next) PetscOptionsObject->next = *amsopt;
100a297a907SKarl Rupp   else {
101e55864a3SBarry Smith     next = PetscOptionsObject->next;
10253acd3b1SBarry Smith     while (next->next) next = next->next;
10353acd3b1SBarry Smith     next->next = *amsopt;
10453acd3b1SBarry Smith   }
10553acd3b1SBarry Smith   PetscFunctionReturn(0);
10653acd3b1SBarry Smith }
10753acd3b1SBarry Smith 
10853acd3b1SBarry Smith #undef __FUNCT__
109aee2cecaSBarry Smith #define __FUNCT__ "PetscScanString"
110aee2cecaSBarry Smith /*
1113fc1eb6aSBarry Smith     PetscScanString -  Gets user input via stdin from process and broadcasts to all processes
1123fc1eb6aSBarry Smith 
1133fc1eb6aSBarry Smith     Collective on MPI_Comm
1143fc1eb6aSBarry Smith 
1153fc1eb6aSBarry Smith    Input Parameters:
1163fc1eb6aSBarry Smith +     commm - communicator for the broadcast, must be PETSC_COMM_WORLD
1173fc1eb6aSBarry Smith .     n - length of the string, must be the same on all processes
1183fc1eb6aSBarry Smith -     str - location to store input
119aee2cecaSBarry Smith 
120aee2cecaSBarry Smith     Bugs:
121aee2cecaSBarry Smith .   Assumes process 0 of the given communicator has access to stdin
122aee2cecaSBarry Smith 
123aee2cecaSBarry Smith */
1243fc1eb6aSBarry Smith static PetscErrorCode PetscScanString(MPI_Comm comm,size_t n,char str[])
125aee2cecaSBarry Smith {
126330cf3c9SBarry Smith   size_t         i;
127aee2cecaSBarry Smith   char           c;
1283fc1eb6aSBarry Smith   PetscMPIInt    rank,nm;
129aee2cecaSBarry Smith   PetscErrorCode ierr;
130aee2cecaSBarry Smith 
131aee2cecaSBarry Smith   PetscFunctionBegin;
132aee2cecaSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
133aee2cecaSBarry Smith   if (!rank) {
134aee2cecaSBarry Smith     c = (char) getchar();
135aee2cecaSBarry Smith     i = 0;
136aee2cecaSBarry Smith     while (c != '\n' && i < n-1) {
137aee2cecaSBarry Smith       str[i++] = c;
138aee2cecaSBarry Smith       c = (char)getchar();
139aee2cecaSBarry Smith     }
140aee2cecaSBarry Smith     str[i] = 0;
141aee2cecaSBarry Smith   }
1424dc2109aSBarry Smith   ierr = PetscMPIIntCast(n,&nm);CHKERRQ(ierr);
1433fc1eb6aSBarry Smith   ierr = MPI_Bcast(str,nm,MPI_CHAR,0,comm);CHKERRQ(ierr);
144aee2cecaSBarry Smith   PetscFunctionReturn(0);
145aee2cecaSBarry Smith }
146aee2cecaSBarry Smith 
147ead66b60SBarry Smith #undef __FUNCT__
148ead66b60SBarry Smith #define __FUNCT__ "PetscStrdup"
1495b02f95dSBarry Smith /*
1505b02f95dSBarry Smith     This is needed because certain strings may be freed by SAWs, hence we cannot use PetscStrallocpy()
1515b02f95dSBarry Smith */
1525b02f95dSBarry Smith static PetscErrorCode  PetscStrdup(const char s[],char *t[])
1535b02f95dSBarry Smith {
1545b02f95dSBarry Smith   PetscErrorCode ierr;
1555b02f95dSBarry Smith   size_t         len;
1565b02f95dSBarry Smith   char           *tmp = 0;
1575b02f95dSBarry Smith 
1585b02f95dSBarry Smith   PetscFunctionBegin;
1595b02f95dSBarry Smith   if (s) {
1605b02f95dSBarry Smith     ierr = PetscStrlen(s,&len);CHKERRQ(ierr);
161ee48884fSBarry Smith     tmp = (char*) malloc((len+1)*sizeof(char*));
1625b02f95dSBarry Smith     if (!tmp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"No memory to duplicate string");
1635b02f95dSBarry Smith     ierr = PetscStrcpy(tmp,s);CHKERRQ(ierr);
1645b02f95dSBarry Smith   }
1655b02f95dSBarry Smith   *t = tmp;
1665b02f95dSBarry Smith   PetscFunctionReturn(0);
1675b02f95dSBarry Smith }
1685b02f95dSBarry Smith 
1695b02f95dSBarry Smith 
170aee2cecaSBarry Smith #undef __FUNCT__
171aee2cecaSBarry Smith #define __FUNCT__ "PetscOptionsGetFromTextInput"
172aee2cecaSBarry Smith /*
1733cc1e11dSBarry Smith     PetscOptionsGetFromTextInput - Presents all the PETSc Options processed by the program so the user may change them at runtime
174aee2cecaSBarry Smith 
175aee2cecaSBarry Smith     Notes: this isn't really practical, it is just to demonstrate the principle
176aee2cecaSBarry Smith 
1777781c08eSBarry Smith     A carriage return indicates no change from the default; but this like -ksp_monitor <stdout>  the default is actually not stdout the default
1787781c08eSBarry Smith     is to do nothing so to get it to use stdout you need to type stdout. This is kind of bug?
1797781c08eSBarry Smith 
180aee2cecaSBarry Smith     Bugs:
1817781c08eSBarry Smith +    All processes must traverse through the exact same set of option queries due to the call to PetscScanString()
1823cc1e11dSBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
183aee2cecaSBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
184aee2cecaSBarry Smith 
1853cc1e11dSBarry Smith     Developer Notes: Normally the GUI that presents the options the user and retrieves the values would be running in a different
1863cc1e11dSBarry Smith      address space and communicating with the PETSc program
1873cc1e11dSBarry Smith 
188aee2cecaSBarry Smith */
1898c34d3f5SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptions *PetscOptionsObject)
1906356e834SBarry Smith {
1916356e834SBarry Smith   PetscErrorCode ierr;
1928c34d3f5SBarry Smith   PetscOption   next = PetscOptionsObject->next;
1936356e834SBarry Smith   char           str[512];
1947781c08eSBarry Smith   PetscBool      bid;
195a4404d99SBarry Smith   PetscReal      ir,*valr;
196330cf3c9SBarry Smith   PetscInt       *vald;
197330cf3c9SBarry Smith   size_t         i;
1986356e834SBarry Smith 
1992a409bb0SBarry Smith   PetscFunctionBegin;
200e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
2016356e834SBarry Smith   while (next) {
2026356e834SBarry Smith     switch (next->type) {
2036356e834SBarry Smith     case OPTION_HEAD:
2046356e834SBarry Smith       break;
205e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
206e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
207e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
208e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
209e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
210e26ddf31SBarry Smith         if (i < next->arraylength-1) {
211e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
212e26ddf31SBarry Smith         }
213e26ddf31SBarry Smith       }
214e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
215e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
216e26ddf31SBarry Smith       if (str[0]) {
217e26ddf31SBarry Smith         PetscToken token;
218e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
219e26ddf31SBarry Smith         size_t     len;
220e26ddf31SBarry Smith         char       *value;
221ace3abfcSBarry Smith         PetscBool  foundrange;
222e26ddf31SBarry Smith 
223e26ddf31SBarry Smith         next->set = PETSC_TRUE;
224e26ddf31SBarry Smith         value     = str;
225e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
226e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
227e26ddf31SBarry Smith         while (n < nmax) {
228e26ddf31SBarry Smith           if (!value) break;
229e26ddf31SBarry Smith 
230e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
231e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
232e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
233e26ddf31SBarry Smith           if (value[0] == '-') i=2;
234e26ddf31SBarry Smith           else i=1;
235330cf3c9SBarry Smith           for (;i<len; i++) {
236e26ddf31SBarry Smith             if (value[i] == '-') {
237e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
238e26ddf31SBarry Smith               value[i] = 0;
239cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
240cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
241e32f2f54SBarry 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);
242e32f2f54SBarry 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);
243e26ddf31SBarry Smith               for (; start<end; start++) {
244e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
245e26ddf31SBarry Smith               }
246e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
247e26ddf31SBarry Smith               break;
248e26ddf31SBarry Smith             }
249e26ddf31SBarry Smith           }
250e26ddf31SBarry Smith           if (!foundrange) {
251cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
252e26ddf31SBarry Smith             dvalue++;
253e26ddf31SBarry Smith             n++;
254e26ddf31SBarry Smith           }
255e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
256e26ddf31SBarry Smith         }
2578c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
258e26ddf31SBarry Smith       }
259e26ddf31SBarry Smith       break;
260e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
261e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
262e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
263e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
264e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
265e26ddf31SBarry Smith         if (i < next->arraylength-1) {
266e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
267e26ddf31SBarry Smith         }
268e26ddf31SBarry Smith       }
269e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
270e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
271e26ddf31SBarry Smith       if (str[0]) {
272e26ddf31SBarry Smith         PetscToken token;
273e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
274e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
275e26ddf31SBarry Smith         char       *value;
276e26ddf31SBarry Smith 
277e26ddf31SBarry Smith         next->set = PETSC_TRUE;
278e26ddf31SBarry Smith         value     = str;
279e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
280e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
281e26ddf31SBarry Smith         while (n < nmax) {
282e26ddf31SBarry Smith           if (!value) break;
283cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
284e26ddf31SBarry Smith           dvalue++;
285e26ddf31SBarry Smith           n++;
286e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
287e26ddf31SBarry Smith         }
2888c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
289e26ddf31SBarry Smith       }
290e26ddf31SBarry Smith       break;
2916356e834SBarry Smith     case OPTION_INT:
292e55864a3SBarry 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);
2933fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2943fc1eb6aSBarry Smith       if (str[0]) {
295d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
296d25d7f95SJed Brown         long long lid;
297d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
2981a1499c8SBarry Smith         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %lld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid);
299c272547aSJed Brown #else
300d25d7f95SJed Brown         long  lid;
301d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
3021a1499c8SBarry Smith         if (lid > PETSC_MAX_INT || lid < PETSC_MIN_INT) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Argument: -%s%s %ld",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,lid);
303c272547aSJed Brown #endif
304a297a907SKarl Rupp 
305d25d7f95SJed Brown         next->set = PETSC_TRUE;
306d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
307aee2cecaSBarry Smith       }
308aee2cecaSBarry Smith       break;
309aee2cecaSBarry Smith     case OPTION_REAL:
310e55864a3SBarry 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);
3113fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3123fc1eb6aSBarry Smith       if (str[0]) {
313ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
314a4404d99SBarry Smith         sscanf(str,"%e",&ir);
315ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
316aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
317ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
318d9822059SBarry Smith         ir = strtoflt128(str,0);
319d9822059SBarry Smith #else
320513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
321a4404d99SBarry Smith #endif
322aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
323aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
324aee2cecaSBarry Smith       }
325aee2cecaSBarry Smith       break;
3267781c08eSBarry Smith     case OPTION_BOOL:
32783355fc5SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <%s>: %s (%s) ",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1,*(PetscBool*)next->data ? "true": "false",next->text,next->man);CHKERRQ(ierr);
3287781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3297781c08eSBarry Smith       if (str[0]) {
3307781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3317781c08eSBarry Smith         next->set = PETSC_TRUE;
3327781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3337781c08eSBarry Smith       }
3347781c08eSBarry Smith       break;
335aee2cecaSBarry Smith     case OPTION_STRING:
336e55864a3SBarry 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);
3373fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3383fc1eb6aSBarry Smith       if (str[0]) {
339aee2cecaSBarry Smith         next->set = PETSC_TRUE;
34064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3415b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3426356e834SBarry Smith       }
3436356e834SBarry Smith       break;
344a264d7a6SBarry Smith     case OPTION_FLIST:
345e55864a3SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3463cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3473cc1e11dSBarry Smith       if (str[0]) {
348e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3493cc1e11dSBarry Smith         next->set = PETSC_TRUE;
35064facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3515b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3523cc1e11dSBarry Smith       }
3533cc1e11dSBarry Smith       break;
354b432afdaSMatthew Knepley     default:
355b432afdaSMatthew Knepley       break;
3566356e834SBarry Smith     }
3576356e834SBarry Smith     next = next->next;
3586356e834SBarry Smith   }
3596356e834SBarry Smith   PetscFunctionReturn(0);
3606356e834SBarry Smith }
3616356e834SBarry Smith 
362e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
363e04113cfSBarry Smith #include <petscviewersaws.h>
364d5649816SBarry Smith 
365d5649816SBarry Smith static int count = 0;
366d5649816SBarry Smith 
367b3506946SBarry Smith #undef __FUNCT__
368e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
369e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
370d5649816SBarry Smith {
3712657e9d9SBarry Smith   PetscFunctionBegin;
372d5649816SBarry Smith   PetscFunctionReturn(0);
373d5649816SBarry Smith }
374d5649816SBarry Smith 
3759c1e0ce8SBarry Smith static const char *OptionsHeader = "<head>\n"
37623a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/jquery-1.9.1.js\"></script>\n"
37723a1ff79SBarry Smith                                    "<script type=\"text/javascript\" src=\"http://www.mcs.anl.gov/research/projects/saws/js/SAWs.js\"></script>\n"
378d1fc0251SBarry Smith                                    "<script type=\"text/javascript\" src=\"js/PETSc.js\"></script>\n"
37964bbc9efSBarry Smith                                    "<script>\n"
38064bbc9efSBarry Smith                                       "jQuery(document).ready(function() {\n"
38164bbc9efSBarry Smith                                          "PETSc.getAndDisplayDirectory(null,\"#variablesInfo\")\n"
38264bbc9efSBarry Smith                                       "})\n"
38364bbc9efSBarry Smith                                   "</script>\n"
38464bbc9efSBarry Smith                                   "</head>\n";
3851423471aSBarry Smith 
3861423471aSBarry Smith /*  Determines the size and style of the scroll region where PETSc options selectable from users are displayed */
3871423471aSBarry Smith static const char *OptionsBodyBottom = "<div id=\"variablesInfo\" style=\"background-color:lightblue;height:auto;max-height:500px;overflow:scroll;\"></div>\n<br>\n</body>";
38864bbc9efSBarry Smith 
389d5649816SBarry Smith #undef __FUNCT__
3907781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
391b3506946SBarry Smith /*
3927781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
393b3506946SBarry Smith 
394b3506946SBarry Smith     Bugs:
395b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
396b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
397b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
398b3506946SBarry Smith 
399b3506946SBarry Smith 
400b3506946SBarry Smith */
401f7b25cf6SBarry Smith PetscErrorCode PetscOptionsSAWsInput(PetscOptions *PetscOptionsObject)
402b3506946SBarry Smith {
403b3506946SBarry Smith   PetscErrorCode ierr;
4048c34d3f5SBarry Smith   PetscOption    next     = PetscOptionsObject->next;
405d5649816SBarry Smith   static int     mancount = 0;
406b3506946SBarry Smith   char           options[16];
4077aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
408a23eb982SSurtai Han   PetscBool      stopasking    = PETSC_FALSE;
40988a9752cSBarry Smith   char           manname[16],textname[16];
4102657e9d9SBarry Smith   char           dir[1024];
411b3506946SBarry Smith 
4122a409bb0SBarry Smith   PetscFunctionBegin;
413b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
414b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
415a297a907SKarl Rupp 
4167325c4a4SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* SAWs will change this, so cannot pass prefix directly */
4171bc75a8dSBarry Smith 
418d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
41983355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4207781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
42183355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4222657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
423a23eb982SSurtai Han   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/StopAsking",&stopasking,1,SAWs_WRITE,SAWs_BOOLEAN));
424b3506946SBarry Smith 
425b3506946SBarry Smith   while (next) {
426d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4272657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4287781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
429d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4307781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4317781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4329f32e415SBarry Smith 
433b3506946SBarry Smith     switch (next->type) {
434b3506946SBarry Smith     case OPTION_HEAD:
435b3506946SBarry Smith       break;
436b3506946SBarry Smith     case OPTION_INT_ARRAY:
4377781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4382657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
439b3506946SBarry Smith       break;
440b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4417781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4422657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
443b3506946SBarry Smith       break;
444b3506946SBarry Smith     case OPTION_INT:
4457781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4462657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
447b3506946SBarry Smith       break;
448b3506946SBarry Smith     case OPTION_REAL:
4497781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4502657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
451b3506946SBarry Smith       break;
4527781c08eSBarry Smith     case OPTION_BOOL:
4537781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4542657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4551ae3d29cSBarry Smith       break;
4567781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4577781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4582657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
45971f08665SBarry Smith       break;
460b3506946SBarry Smith     case OPTION_STRING:
4617781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4627781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4631ae3d29cSBarry Smith       break;
4641ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4657781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4662657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
467b3506946SBarry Smith       break;
468a264d7a6SBarry Smith     case OPTION_FLIST:
469a264d7a6SBarry Smith       {
470a264d7a6SBarry Smith       PetscInt ntext;
4717781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4727781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
473a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
474a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
475a264d7a6SBarry Smith       }
476a264d7a6SBarry Smith       break;
4771ae3d29cSBarry Smith     case OPTION_ELIST:
478a264d7a6SBarry Smith       {
479a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4807781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4817781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
482ead66b60SBarry Smith       ierr = PetscMalloc1((ntext+1),(char***)&next->edata);CHKERRQ(ierr);
4831ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
484a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
485a264d7a6SBarry Smith       }
486a264d7a6SBarry Smith       break;
487b3506946SBarry Smith     default:
488b3506946SBarry Smith       break;
489b3506946SBarry Smith     }
490b3506946SBarry Smith     next = next->next;
491b3506946SBarry Smith   }
492b3506946SBarry Smith 
493b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
49464bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Header,("index.html",OptionsHeader));
49564bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Push_Body,("index.html",2,OptionsBodyBottom));
4967aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
49764bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Header,("index.html"));
49864bbc9efSBarry Smith   PetscStackCallSAWs(SAWs_Pop_Body,("index.html",2));
499b3506946SBarry Smith 
50088a9752cSBarry Smith   /* determine if any values have been set in GUI */
50183355fc5SBarry Smith   next = PetscOptionsObject->next;
50288a9752cSBarry Smith   while (next) {
50388a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
504f7b25cf6SBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,(int*)&next->set));
50588a9752cSBarry Smith     next = next->next;
50688a9752cSBarry Smith   }
50788a9752cSBarry Smith 
508b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
509f7b25cf6SBarry Smith   if (changedmethod) PetscOptionsObject->count = -2;
510b3506946SBarry Smith 
511a23eb982SSurtai Han   if (stopasking) {
512a23eb982SSurtai Han     PetscOptionsPublish      = PETSC_FALSE;
51312655325SBarry Smith     PetscOptionsObject->count = 0;//do not ask for same thing again
514a23eb982SSurtai Han   }
515a23eb982SSurtai Han 
5169a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
517b3506946SBarry Smith   PetscFunctionReturn(0);
518b3506946SBarry Smith }
519b3506946SBarry Smith #endif
520b3506946SBarry Smith 
5216356e834SBarry Smith #undef __FUNCT__
52253acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
5238c34d3f5SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptions *PetscOptionsObject)
52453acd3b1SBarry Smith {
52553acd3b1SBarry Smith   PetscErrorCode ierr;
5268c34d3f5SBarry Smith   PetscOption   last;
5276356e834SBarry Smith   char           option[256],value[1024],tmp[32];
528330cf3c9SBarry Smith   size_t         j;
52953acd3b1SBarry Smith 
53053acd3b1SBarry Smith   PetscFunctionBegin;
53183355fc5SBarry Smith   if (PetscOptionsObject->next) {
53283355fc5SBarry Smith     if (!PetscOptionsObject->count) {
533a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
534f7b25cf6SBarry Smith       ierr = PetscOptionsSAWsInput(PetscOptionsObject);CHKERRQ(ierr);
535b3506946SBarry Smith #else
536e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
537b3506946SBarry Smith #endif
538aee2cecaSBarry Smith     }
539aee2cecaSBarry Smith   }
5406356e834SBarry Smith 
541e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
5426356e834SBarry Smith 
543e26ddf31SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
544e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
5457a72a596SBarry Smith   /* reset alreadyprinted flag */
546e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
547e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
548e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
54953acd3b1SBarry Smith 
550e55864a3SBarry Smith   while (PetscOptionsObject->next) {
551e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
552e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
55353acd3b1SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
554e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
555e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5566356e834SBarry Smith       } else {
557e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5586356e834SBarry Smith       }
5596356e834SBarry Smith 
560e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5616356e834SBarry Smith       case OPTION_HEAD:
5626356e834SBarry Smith         break;
5636356e834SBarry Smith       case OPTION_INT_ARRAY:
564e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
565e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
566e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
5676356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5686356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5696356e834SBarry Smith         }
5706356e834SBarry Smith         break;
5716356e834SBarry Smith       case OPTION_INT:
572e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5736356e834SBarry Smith         break;
5746356e834SBarry Smith       case OPTION_REAL:
575e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5766356e834SBarry Smith         break;
5776356e834SBarry Smith       case OPTION_REAL_ARRAY:
578e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
579e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
580e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5816356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5826356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5836356e834SBarry Smith         }
5846356e834SBarry Smith         break;
585*050cccc3SHong Zhang       case OPTION_SCALAR_ARRAY:
586*050cccc3SHong Zhang         sprintf(value,"%g+%gi",(double)PetscRealPart(((PetscReal*)PetscOptionsObject->next->data)[0]),(double)PetscImaginaryPart(((PetscReal*)PetscOptionsObject->next->data)[0]));
587*050cccc3SHong Zhang         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
588*050cccc3SHong Zhang           sprintf(tmp,"%g+%gi",(double)PetscRealPart(((PetscReal*)PetscOptionsObject->next->data)[j]),(double)PetscImaginaryPart(((PetscReal*)PetscOptionsObject->next->data)[j]));
589*050cccc3SHong Zhang           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
590*050cccc3SHong Zhang           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
591*050cccc3SHong Zhang         }
592*050cccc3SHong Zhang         break;
5937781c08eSBarry Smith       case OPTION_BOOL:
594e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5956356e834SBarry Smith         break;
5967781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
597e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
598e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
599e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
6001ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6011ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6021ae3d29cSBarry Smith         }
6031ae3d29cSBarry Smith         break;
604a264d7a6SBarry Smith       case OPTION_FLIST:
6056991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6066991f827SBarry Smith         break;
6071ae3d29cSBarry Smith       case OPTION_ELIST:
608e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6096356e834SBarry Smith         break;
6101ae3d29cSBarry Smith       case OPTION_STRING:
611e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
612d51da6bfSBarry Smith         break;
6131ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
614e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
615e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
616e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6171ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6181ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6191ae3d29cSBarry Smith         }
6206356e834SBarry Smith         break;
6216356e834SBarry Smith       }
6226356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6236356e834SBarry Smith     }
6246991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
6256991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
6266991f827SBarry Smith     }
627e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
628e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
629e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
630e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
631c979a496SBarry Smith 
63283355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
63383355fc5SBarry Smith       free(PetscOptionsObject->next->data);
634c979a496SBarry Smith     } else {
63583355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
636c979a496SBarry Smith     }
6377781c08eSBarry Smith 
63883355fc5SBarry Smith     last                    = PetscOptionsObject->next;
63983355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6406356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6416356e834SBarry Smith   }
642f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
643e55864a3SBarry Smith   PetscOptionsObject->next = 0;
64453acd3b1SBarry Smith   PetscFunctionReturn(0);
64553acd3b1SBarry Smith }
64653acd3b1SBarry Smith 
64753acd3b1SBarry Smith #undef __FUNCT__
648e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
64953acd3b1SBarry Smith /*@C
65053acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
65153acd3b1SBarry Smith 
6523f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
65353acd3b1SBarry Smith 
65453acd3b1SBarry Smith    Input Parameters:
65553acd3b1SBarry Smith +  opt - option name
65653acd3b1SBarry Smith .  text - short string that describes the option
65753acd3b1SBarry Smith .  man - manual page with additional information on option
65853acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6590fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6600fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6610fdccdaeSBarry Smith $                 value = defaultvalue
6620fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6630fdccdaeSBarry Smith $                 if (flg) {
66453acd3b1SBarry Smith 
66553acd3b1SBarry Smith    Output Parameter:
66653acd3b1SBarry Smith +  value - the  value to return
667b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
66853acd3b1SBarry Smith 
66953acd3b1SBarry Smith    Level: beginner
67053acd3b1SBarry Smith 
67153acd3b1SBarry Smith    Concepts: options database
67253acd3b1SBarry Smith 
67353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
67453acd3b1SBarry Smith 
67553acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
67653acd3b1SBarry Smith 
67753acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
678acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
679acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
68053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
68153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
682acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
683a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
68453acd3b1SBarry Smith @*/
6858c34d3f5SBarry Smith PetscErrorCode  PetscOptionsEnum_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
68653acd3b1SBarry Smith {
68753acd3b1SBarry Smith   PetscErrorCode ierr;
68853acd3b1SBarry Smith   PetscInt       ntext = 0;
689aa5bb8c0SSatish Balay   PetscInt       tval;
690ace3abfcSBarry Smith   PetscBool      tflg;
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith   PetscFunctionBegin;
69353acd3b1SBarry Smith   while (list[ntext++]) {
694e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
69553acd3b1SBarry Smith   }
696e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
69753acd3b1SBarry Smith   ntext -= 3;
698e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
699aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
700aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
701aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
70253acd3b1SBarry Smith   PetscFunctionReturn(0);
70353acd3b1SBarry Smith }
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
70653acd3b1SBarry Smith #undef __FUNCT__
707e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
70853acd3b1SBarry Smith /*@C
70953acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
71053acd3b1SBarry Smith 
7113f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
71253acd3b1SBarry Smith 
71353acd3b1SBarry Smith    Input Parameters:
71453acd3b1SBarry Smith +  opt - option name
71553acd3b1SBarry Smith .  text - short string that describes the option
71653acd3b1SBarry Smith .  man - manual page with additional information on option
7170fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7180fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7190fdccdaeSBarry Smith $                 value = defaultvalue
7200fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7210fdccdaeSBarry Smith $                 if (flg) {
72253acd3b1SBarry Smith 
72353acd3b1SBarry Smith    Output Parameter:
72453acd3b1SBarry Smith +  value - the integer value to return
72553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
72653acd3b1SBarry Smith 
72753acd3b1SBarry Smith    Level: beginner
72853acd3b1SBarry Smith 
72953acd3b1SBarry Smith    Concepts: options database^has int
73053acd3b1SBarry Smith 
73153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
73253acd3b1SBarry Smith 
73353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
734acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
735acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
73653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
73753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
738acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
739a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
74053acd3b1SBarry Smith @*/
7418c34d3f5SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
74253acd3b1SBarry Smith {
74353acd3b1SBarry Smith   PetscErrorCode ierr;
7448c34d3f5SBarry Smith   PetscOption    amsopt;
74512655325SBarry Smith   PetscBool      wasset;
74653acd3b1SBarry Smith 
74753acd3b1SBarry Smith   PetscFunctionBegin;
748e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
749e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7506356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
75112655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
7523e211508SBarry Smith 
75312655325SBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
7543e211508SBarry Smith     if (wasset) {
75512655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
7563e211508SBarry Smith     }
757af6d86caSBarry Smith   }
758e55864a3SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
759e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
7601a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
76153acd3b1SBarry Smith   }
76253acd3b1SBarry Smith   PetscFunctionReturn(0);
76353acd3b1SBarry Smith }
76453acd3b1SBarry Smith 
76553acd3b1SBarry Smith #undef __FUNCT__
766e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
76753acd3b1SBarry Smith /*@C
76853acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
76953acd3b1SBarry Smith 
7703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
77153acd3b1SBarry Smith 
77253acd3b1SBarry Smith    Input Parameters:
77353acd3b1SBarry Smith +  opt - option name
77453acd3b1SBarry Smith .  text - short string that describes the option
77553acd3b1SBarry Smith .  man - manual page with additional information on option
7760fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
777bcbf2dc5SJed Brown -  len - length of the result string including null terminator
77853acd3b1SBarry Smith 
77953acd3b1SBarry Smith    Output Parameter:
78053acd3b1SBarry Smith +  value - the value to return
78153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
78253acd3b1SBarry Smith 
78353acd3b1SBarry Smith    Level: beginner
78453acd3b1SBarry Smith 
78553acd3b1SBarry Smith    Concepts: options database^has int
78653acd3b1SBarry Smith 
78753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
78853acd3b1SBarry Smith 
7897fccdfe4SBarry 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).
7907fccdfe4SBarry Smith 
79153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
792acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
793acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
79453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
79553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
796acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
797a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
79853acd3b1SBarry Smith @*/
7998c34d3f5SBarry Smith PetscErrorCode  PetscOptionsString_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
80053acd3b1SBarry Smith {
80153acd3b1SBarry Smith   PetscErrorCode ierr;
8028c34d3f5SBarry Smith   PetscOption   amsopt;
80353acd3b1SBarry Smith 
80453acd3b1SBarry Smith   PetscFunctionBegin;
8051a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
8061a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
80764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8080fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
809af6d86caSBarry Smith   }
810e55864a3SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
811e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8121a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
81353acd3b1SBarry Smith   }
81453acd3b1SBarry Smith   PetscFunctionReturn(0);
81553acd3b1SBarry Smith }
81653acd3b1SBarry Smith 
81753acd3b1SBarry Smith #undef __FUNCT__
818e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
81953acd3b1SBarry Smith /*@C
82053acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
82153acd3b1SBarry Smith 
8223f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
82353acd3b1SBarry Smith 
82453acd3b1SBarry Smith    Input Parameters:
82553acd3b1SBarry Smith +  opt - option name
82653acd3b1SBarry Smith .  text - short string that describes the option
82753acd3b1SBarry Smith .  man - manual page with additional information on option
8280fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8290fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
8300fdccdaeSBarry Smith $                 value = defaultvalue
8310fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
8320fdccdaeSBarry Smith $                 if (flg) {
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Output Parameter:
83553acd3b1SBarry Smith +  value - the value to return
83653acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
83753acd3b1SBarry Smith 
83853acd3b1SBarry Smith    Level: beginner
83953acd3b1SBarry Smith 
84053acd3b1SBarry Smith    Concepts: options database^has int
84153acd3b1SBarry Smith 
84253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
84353acd3b1SBarry Smith 
84453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
845acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
846acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
84753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
849acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
850a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
85153acd3b1SBarry Smith @*/
8528c34d3f5SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
85353acd3b1SBarry Smith {
85453acd3b1SBarry Smith   PetscErrorCode ierr;
8558c34d3f5SBarry Smith   PetscOption   amsopt;
85653acd3b1SBarry Smith 
85753acd3b1SBarry Smith   PetscFunctionBegin;
858e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
859e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
860538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
861a297a907SKarl Rupp 
8620fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
863538aa990SBarry Smith   }
8641a1499c8SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
8651a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8661a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
86753acd3b1SBarry Smith   }
86853acd3b1SBarry Smith   PetscFunctionReturn(0);
86953acd3b1SBarry Smith }
87053acd3b1SBarry Smith 
87153acd3b1SBarry Smith #undef __FUNCT__
872e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
87353acd3b1SBarry Smith /*@C
87453acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
87553acd3b1SBarry Smith 
8763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
87753acd3b1SBarry Smith 
87853acd3b1SBarry Smith    Input Parameters:
87953acd3b1SBarry Smith +  opt - option name
88053acd3b1SBarry Smith .  text - short string that describes the option
88153acd3b1SBarry Smith .  man - manual page with additional information on option
8820fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8830fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
8840fdccdaeSBarry Smith $                 value = defaultvalue
8850fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
8860fdccdaeSBarry Smith $                 if (flg) {
8870fdccdaeSBarry Smith 
88853acd3b1SBarry Smith 
88953acd3b1SBarry Smith    Output Parameter:
89053acd3b1SBarry Smith +  value - the value to return
89153acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
89253acd3b1SBarry Smith 
89353acd3b1SBarry Smith    Level: beginner
89453acd3b1SBarry Smith 
89553acd3b1SBarry Smith    Concepts: options database^has int
89653acd3b1SBarry Smith 
89753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89853acd3b1SBarry Smith 
89953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
900acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
901acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
90253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
904acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
905a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
90653acd3b1SBarry Smith @*/
9078c34d3f5SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
90853acd3b1SBarry Smith {
90953acd3b1SBarry Smith   PetscErrorCode ierr;
91053acd3b1SBarry Smith 
91153acd3b1SBarry Smith   PetscFunctionBegin;
91253acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9130fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
91453acd3b1SBarry Smith #else
915e55864a3SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
91653acd3b1SBarry Smith #endif
91753acd3b1SBarry Smith   PetscFunctionReturn(0);
91853acd3b1SBarry Smith }
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith #undef __FUNCT__
921e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
92253acd3b1SBarry Smith /*@C
92390d69ab7SBarry 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
92490d69ab7SBarry Smith                       its value is set to false.
92553acd3b1SBarry Smith 
9263f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
92753acd3b1SBarry Smith 
92853acd3b1SBarry Smith    Input Parameters:
92953acd3b1SBarry Smith +  opt - option name
93053acd3b1SBarry Smith .  text - short string that describes the option
93153acd3b1SBarry Smith -  man - manual page with additional information on option
93253acd3b1SBarry Smith 
93353acd3b1SBarry Smith    Output Parameter:
93453acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
93553acd3b1SBarry Smith 
93653acd3b1SBarry Smith    Level: beginner
93753acd3b1SBarry Smith 
93853acd3b1SBarry Smith    Concepts: options database^has int
93953acd3b1SBarry Smith 
94053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
94153acd3b1SBarry Smith 
94253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
943acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
944acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
94553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
94653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
947acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
948a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
94953acd3b1SBarry Smith @*/
9508c34d3f5SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
95153acd3b1SBarry Smith {
95253acd3b1SBarry Smith   PetscErrorCode ierr;
9538c34d3f5SBarry Smith   PetscOption   amsopt;
95453acd3b1SBarry Smith 
95553acd3b1SBarry Smith   PetscFunctionBegin;
956e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
957e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
958ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
959a297a907SKarl Rupp 
960ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9611ae3d29cSBarry Smith   }
962e55864a3SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
963e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
964e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
96553acd3b1SBarry Smith   }
96653acd3b1SBarry Smith   PetscFunctionReturn(0);
96753acd3b1SBarry Smith }
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith #undef __FUNCT__
970e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
97153acd3b1SBarry Smith /*@C
972a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
97353acd3b1SBarry Smith 
9743f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
97553acd3b1SBarry Smith 
97653acd3b1SBarry Smith    Input Parameters:
97753acd3b1SBarry Smith +  opt - option name
97853acd3b1SBarry Smith .  text - short string that describes the option
97953acd3b1SBarry Smith .  man - manual page with additional information on option
98053acd3b1SBarry Smith .  list - the possible choices
9810fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
9820fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
9830fdccdaeSBarry Smith $                 if (flg) {
9843cc1e11dSBarry Smith -  len - the length of the character array value
98553acd3b1SBarry Smith 
98653acd3b1SBarry Smith    Output Parameter:
98753acd3b1SBarry Smith +  value - the value to return
98853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
98953acd3b1SBarry Smith 
99053acd3b1SBarry Smith    Level: intermediate
99153acd3b1SBarry Smith 
99253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
99353acd3b1SBarry Smith 
99453acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
99553acd3b1SBarry Smith 
99653acd3b1SBarry Smith    To get a listing of all currently specified options,
99788c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
99853acd3b1SBarry Smith 
999eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
1000eabe10d7SBarry Smith 
100153acd3b1SBarry Smith    Concepts: options database^list
100253acd3b1SBarry Smith 
100353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1004acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
100553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
100653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1007acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1008a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
100953acd3b1SBarry Smith @*/
10108c34d3f5SBarry Smith PetscErrorCode  PetscOptionsFList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
101153acd3b1SBarry Smith {
101253acd3b1SBarry Smith   PetscErrorCode ierr;
10138c34d3f5SBarry Smith   PetscOption   amsopt;
101453acd3b1SBarry Smith 
101553acd3b1SBarry Smith   PetscFunctionBegin;
10161a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10171a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
101864facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10190fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10203cc1e11dSBarry Smith     amsopt->flist = list;
10213cc1e11dSBarry Smith   }
10221a1499c8SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
10231a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10241a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
102553acd3b1SBarry Smith   }
102653acd3b1SBarry Smith   PetscFunctionReturn(0);
102753acd3b1SBarry Smith }
102853acd3b1SBarry Smith 
102953acd3b1SBarry Smith #undef __FUNCT__
1030e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
103153acd3b1SBarry Smith /*@C
103253acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
103353acd3b1SBarry Smith 
10343f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
103553acd3b1SBarry Smith 
103653acd3b1SBarry Smith    Input Parameters:
103753acd3b1SBarry Smith +  opt - option name
103853acd3b1SBarry Smith .  ltext - short string that describes the option
103953acd3b1SBarry Smith .  man - manual page with additional information on option
1040a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
104153acd3b1SBarry Smith .  ntext - number of choices
10420fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10430fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
10440fdccdaeSBarry Smith $                 if (flg) {
10450fdccdaeSBarry Smith 
104653acd3b1SBarry Smith 
104753acd3b1SBarry Smith    Output Parameter:
104853acd3b1SBarry Smith +  value - the index of the value to return
104953acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
105053acd3b1SBarry Smith 
105153acd3b1SBarry Smith    Level: intermediate
105253acd3b1SBarry Smith 
105353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
105453acd3b1SBarry Smith 
1055a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
105653acd3b1SBarry Smith 
105753acd3b1SBarry Smith    Concepts: options database^list
105853acd3b1SBarry Smith 
105953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1060acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
106153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
106253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1063acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1064a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
106553acd3b1SBarry Smith @*/
10668c34d3f5SBarry Smith PetscErrorCode  PetscOptionsEList_Private(PetscOptions *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
106753acd3b1SBarry Smith {
106853acd3b1SBarry Smith   PetscErrorCode ierr;
106953acd3b1SBarry Smith   PetscInt       i;
10708c34d3f5SBarry Smith   PetscOption   amsopt;
107153acd3b1SBarry Smith 
107253acd3b1SBarry Smith   PetscFunctionBegin;
10731a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10741a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
107564facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10760fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10776991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
10781ae3d29cSBarry Smith     amsopt->nlist = ntext;
10791ae3d29cSBarry Smith   }
10801a1499c8SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
10811a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10821a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
108353acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1084e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
108553acd3b1SBarry Smith     }
1086e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
108753acd3b1SBarry Smith   }
108853acd3b1SBarry Smith   PetscFunctionReturn(0);
108953acd3b1SBarry Smith }
109053acd3b1SBarry Smith 
109153acd3b1SBarry Smith #undef __FUNCT__
1092e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
109353acd3b1SBarry Smith /*@C
1094acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1095d5649816SBarry Smith        which at most a single value can be true.
109653acd3b1SBarry Smith 
10973f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
109853acd3b1SBarry Smith 
109953acd3b1SBarry Smith    Input Parameters:
110053acd3b1SBarry Smith +  opt - option name
110153acd3b1SBarry Smith .  text - short string that describes the option
110253acd3b1SBarry Smith -  man - manual page with additional information on option
110353acd3b1SBarry Smith 
110453acd3b1SBarry Smith    Output Parameter:
110553acd3b1SBarry Smith .  flg - whether that option was set or not
110653acd3b1SBarry Smith 
110753acd3b1SBarry Smith    Level: intermediate
110853acd3b1SBarry Smith 
110953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111053acd3b1SBarry Smith 
1111acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith     Concepts: options database^logical group
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1116acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
111753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
111853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1119acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1120a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
112153acd3b1SBarry Smith @*/
11228c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
112353acd3b1SBarry Smith {
112453acd3b1SBarry Smith   PetscErrorCode ierr;
11258c34d3f5SBarry Smith   PetscOption   amsopt;
112653acd3b1SBarry Smith 
112753acd3b1SBarry Smith   PetscFunctionBegin;
1128e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
112983355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1130ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1131a297a907SKarl Rupp 
1132ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11331ae3d29cSBarry Smith   }
113468b16fdaSBarry Smith   *flg = PETSC_FALSE;
1135e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1136e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1137e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1138e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
113953acd3b1SBarry Smith   }
114053acd3b1SBarry Smith   PetscFunctionReturn(0);
114153acd3b1SBarry Smith }
114253acd3b1SBarry Smith 
114353acd3b1SBarry Smith #undef __FUNCT__
1144e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
114553acd3b1SBarry Smith /*@C
1146acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1147d5649816SBarry Smith        which at most a single value can be true.
114853acd3b1SBarry Smith 
11493f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115053acd3b1SBarry Smith 
115153acd3b1SBarry Smith    Input Parameters:
115253acd3b1SBarry Smith +  opt - option name
115353acd3b1SBarry Smith .  text - short string that describes the option
115453acd3b1SBarry Smith -  man - manual page with additional information on option
115553acd3b1SBarry Smith 
115653acd3b1SBarry Smith    Output Parameter:
115753acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith    Level: intermediate
116053acd3b1SBarry Smith 
116153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116253acd3b1SBarry Smith 
1163acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
116453acd3b1SBarry Smith 
116553acd3b1SBarry Smith     Concepts: options database^logical group
116653acd3b1SBarry Smith 
116753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1168acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
116953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1171acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1172a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
117353acd3b1SBarry Smith @*/
11748c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
117553acd3b1SBarry Smith {
117653acd3b1SBarry Smith   PetscErrorCode ierr;
11778c34d3f5SBarry Smith   PetscOption   amsopt;
117853acd3b1SBarry Smith 
117953acd3b1SBarry Smith   PetscFunctionBegin;
1180e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
118183355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1182ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1183a297a907SKarl Rupp 
1184ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11851ae3d29cSBarry Smith   }
118617326d04SJed Brown   *flg = PETSC_FALSE;
1187e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1188e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1189e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
119053acd3b1SBarry Smith   }
119153acd3b1SBarry Smith   PetscFunctionReturn(0);
119253acd3b1SBarry Smith }
119353acd3b1SBarry Smith 
119453acd3b1SBarry Smith #undef __FUNCT__
1195e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
119653acd3b1SBarry Smith /*@C
1197acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1198d5649816SBarry Smith        which at most a single value can be true.
119953acd3b1SBarry Smith 
12003f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120153acd3b1SBarry Smith 
120253acd3b1SBarry Smith    Input Parameters:
120353acd3b1SBarry Smith +  opt - option name
120453acd3b1SBarry Smith .  text - short string that describes the option
120553acd3b1SBarry Smith -  man - manual page with additional information on option
120653acd3b1SBarry Smith 
120753acd3b1SBarry Smith    Output Parameter:
120853acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
120953acd3b1SBarry Smith 
121053acd3b1SBarry Smith    Level: intermediate
121153acd3b1SBarry Smith 
121253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
121353acd3b1SBarry Smith 
1214acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
121553acd3b1SBarry Smith 
121653acd3b1SBarry Smith     Concepts: options database^logical group
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1219acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
122053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1222acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1223a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
122453acd3b1SBarry Smith @*/
12258c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
122653acd3b1SBarry Smith {
122753acd3b1SBarry Smith   PetscErrorCode ierr;
12288c34d3f5SBarry Smith   PetscOption   amsopt;
122953acd3b1SBarry Smith 
123053acd3b1SBarry Smith   PetscFunctionBegin;
1231e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
123283355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1233ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1234a297a907SKarl Rupp 
1235ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12361ae3d29cSBarry Smith   }
123717326d04SJed Brown   *flg = PETSC_FALSE;
1238e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1239e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1240e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
124153acd3b1SBarry Smith   }
124253acd3b1SBarry Smith   PetscFunctionReturn(0);
124353acd3b1SBarry Smith }
124453acd3b1SBarry Smith 
124553acd3b1SBarry Smith #undef __FUNCT__
1246e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
124753acd3b1SBarry Smith /*@C
1248acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
124953acd3b1SBarry Smith 
12503f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
125153acd3b1SBarry Smith 
125253acd3b1SBarry Smith    Input Parameters:
125353acd3b1SBarry Smith +  opt - option name
125453acd3b1SBarry Smith .  text - short string that describes the option
1255868c398cSBarry Smith .  man - manual page with additional information on option
125694ae4db5SBarry Smith -  currentvalue - the current value
125753acd3b1SBarry Smith 
125853acd3b1SBarry Smith    Output Parameter:
125953acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
126053acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
126153acd3b1SBarry Smith 
126253acd3b1SBarry Smith    Level: beginner
126353acd3b1SBarry Smith 
126453acd3b1SBarry Smith    Concepts: options database^logical
126553acd3b1SBarry Smith 
126653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
126753acd3b1SBarry Smith 
126853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1269acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1270acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
127153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
127253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1273acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1274a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
127553acd3b1SBarry Smith @*/
12768c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
127753acd3b1SBarry Smith {
127853acd3b1SBarry Smith   PetscErrorCode ierr;
1279ace3abfcSBarry Smith   PetscBool      iset;
12808c34d3f5SBarry Smith   PetscOption   amsopt;
128153acd3b1SBarry Smith 
128253acd3b1SBarry Smith   PetscFunctionBegin;
1283e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1284e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1285ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1286a297a907SKarl Rupp 
128794ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1288af6d86caSBarry Smith   }
12891a1499c8SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
129053acd3b1SBarry Smith   if (set) *set = iset;
12911a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
129294ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
12931a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
129453acd3b1SBarry Smith   }
129553acd3b1SBarry Smith   PetscFunctionReturn(0);
129653acd3b1SBarry Smith }
129753acd3b1SBarry Smith 
129853acd3b1SBarry Smith #undef __FUNCT__
1299e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
130053acd3b1SBarry Smith /*@C
130153acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
130253acd3b1SBarry Smith    option in the database. The values must be separated with commas with
130353acd3b1SBarry Smith    no intervening spaces.
130453acd3b1SBarry Smith 
13053f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
130653acd3b1SBarry Smith 
130753acd3b1SBarry Smith    Input Parameters:
130853acd3b1SBarry Smith +  opt - the option one is seeking
130953acd3b1SBarry Smith .  text - short string describing option
131053acd3b1SBarry Smith .  man - manual page for option
131153acd3b1SBarry Smith -  nmax - maximum number of values
131253acd3b1SBarry Smith 
131353acd3b1SBarry Smith    Output Parameter:
131453acd3b1SBarry Smith +  value - location to copy values
131553acd3b1SBarry Smith .  nmax - actual number of values found
131653acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
131753acd3b1SBarry Smith 
131853acd3b1SBarry Smith    Level: beginner
131953acd3b1SBarry Smith 
132053acd3b1SBarry Smith    Notes:
132153acd3b1SBarry Smith    The user should pass in an array of doubles
132253acd3b1SBarry Smith 
132353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
132453acd3b1SBarry Smith 
132553acd3b1SBarry Smith    Concepts: options database^array of strings
132653acd3b1SBarry Smith 
132753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1328acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
132953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
133053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1331acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1332a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
133353acd3b1SBarry Smith @*/
13348c34d3f5SBarry Smith PetscErrorCode PetscOptionsRealArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
133553acd3b1SBarry Smith {
133653acd3b1SBarry Smith   PetscErrorCode ierr;
133753acd3b1SBarry Smith   PetscInt       i;
13388c34d3f5SBarry Smith   PetscOption   amsopt;
133953acd3b1SBarry Smith 
134053acd3b1SBarry Smith   PetscFunctionBegin;
1341e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1342e26ddf31SBarry Smith     PetscReal *vals;
1343e26ddf31SBarry Smith 
1344e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1345e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1346e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1347e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1348e26ddf31SBarry Smith     amsopt->arraylength = *n;
1349e26ddf31SBarry Smith   }
1350e55864a3SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1351e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1352a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
135353acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1354e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
135553acd3b1SBarry Smith     }
1356e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
135753acd3b1SBarry Smith   }
135853acd3b1SBarry Smith   PetscFunctionReturn(0);
135953acd3b1SBarry Smith }
136053acd3b1SBarry Smith 
1361*050cccc3SHong Zhang #undef __FUNCT__
1362*050cccc3SHong Zhang #define __FUNCT__ "PetscOptionsScalarArray_Private"
1363*050cccc3SHong Zhang /*@C
1364*050cccc3SHong Zhang    PetscOptionsScalarArray - Gets an array of Scalar values for a particular
1365*050cccc3SHong Zhang    option in the database. The values must be separated with commas with
1366*050cccc3SHong Zhang    no intervening spaces.
1367*050cccc3SHong Zhang 
1368*050cccc3SHong Zhang    Logically Collective on the communicator passed in PetscOptionsBegin()
1369*050cccc3SHong Zhang 
1370*050cccc3SHong Zhang    Input Parameters:
1371*050cccc3SHong Zhang +  opt - the option one is seeking
1372*050cccc3SHong Zhang .  text - short string describing option
1373*050cccc3SHong Zhang .  man - manual page for option
1374*050cccc3SHong Zhang -  nmax - maximum number of values
1375*050cccc3SHong Zhang 
1376*050cccc3SHong Zhang    Output Parameter:
1377*050cccc3SHong Zhang +  value - location to copy values
1378*050cccc3SHong Zhang .  nmax - actual number of values found
1379*050cccc3SHong Zhang -  set - PETSC_TRUE if found, else PETSC_FALSE
1380*050cccc3SHong Zhang 
1381*050cccc3SHong Zhang    Level: beginner
1382*050cccc3SHong Zhang 
1383*050cccc3SHong Zhang    Notes:
1384*050cccc3SHong Zhang    The user should pass in an array of doubles
1385*050cccc3SHong Zhang 
1386*050cccc3SHong Zhang    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1387*050cccc3SHong Zhang 
1388*050cccc3SHong Zhang    Concepts: options database^array of strings
1389*050cccc3SHong Zhang 
1390*050cccc3SHong Zhang .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1391*050cccc3SHong Zhang           PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1392*050cccc3SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1393*050cccc3SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1394*050cccc3SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1395*050cccc3SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1396*050cccc3SHong Zhang @*/
1397*050cccc3SHong Zhang PetscErrorCode PetscOptionsScalarArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar value[],PetscInt *n,PetscBool  *set)
1398*050cccc3SHong Zhang {
1399*050cccc3SHong Zhang   PetscErrorCode ierr;
1400*050cccc3SHong Zhang   PetscInt       i;
1401*050cccc3SHong Zhang   PetscOption   amsopt;
1402*050cccc3SHong Zhang 
1403*050cccc3SHong Zhang   PetscFunctionBegin;
1404*050cccc3SHong Zhang   if (!PetscOptionsObject->count) {
1405*050cccc3SHong Zhang     PetscScalar *vals;
1406*050cccc3SHong Zhang 
1407*050cccc3SHong Zhang     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_SCALAR_ARRAY,&amsopt);CHKERRQ(ierr);
1408*050cccc3SHong Zhang     ierr = PetscMalloc((*n)*sizeof(PetscScalar),&amsopt->data);CHKERRQ(ierr);
1409*050cccc3SHong Zhang     vals = (PetscScalar*)amsopt->data;
1410*050cccc3SHong Zhang     for (i=0; i<*n; i++) vals[i] = value[i];
1411*050cccc3SHong Zhang     amsopt->arraylength = *n;
1412*050cccc3SHong Zhang   }
1413*050cccc3SHong Zhang   ierr = PetscOptionsGetScalarArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1414*050cccc3SHong Zhang   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1415*050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
1416*050cccc3SHong Zhang     for (i=1; i<*n; i++) {
1417*050cccc3SHong Zhang       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
1418*050cccc3SHong Zhang     }
1419*050cccc3SHong Zhang     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1420*050cccc3SHong Zhang   }
1421*050cccc3SHong Zhang   PetscFunctionReturn(0);
1422*050cccc3SHong Zhang }
142353acd3b1SBarry Smith 
142453acd3b1SBarry Smith #undef __FUNCT__
1425e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
142653acd3b1SBarry Smith /*@C
142753acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1428b32a342fSShri Abhyankar    option in the database.
142953acd3b1SBarry Smith 
14303f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
143153acd3b1SBarry Smith 
143253acd3b1SBarry Smith    Input Parameters:
143353acd3b1SBarry Smith +  opt - the option one is seeking
143453acd3b1SBarry Smith .  text - short string describing option
143553acd3b1SBarry Smith .  man - manual page for option
1436f8a50e2bSBarry Smith -  n - maximum number of values
143753acd3b1SBarry Smith 
143853acd3b1SBarry Smith    Output Parameter:
143953acd3b1SBarry Smith +  value - location to copy values
1440f8a50e2bSBarry Smith .  n - actual number of values found
144153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
144253acd3b1SBarry Smith 
144353acd3b1SBarry Smith    Level: beginner
144453acd3b1SBarry Smith 
144553acd3b1SBarry Smith    Notes:
1446b32a342fSShri Abhyankar    The array can be passed as
1447b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
14480fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
14490fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
14500fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1451b32a342fSShri Abhyankar 
1452b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
145353acd3b1SBarry Smith 
145453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
145553acd3b1SBarry Smith 
1456b32a342fSShri Abhyankar    Concepts: options database^array of ints
145753acd3b1SBarry Smith 
145853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1459acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
146053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
146153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1462acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1463a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
146453acd3b1SBarry Smith @*/
14658c34d3f5SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
146653acd3b1SBarry Smith {
146753acd3b1SBarry Smith   PetscErrorCode ierr;
146853acd3b1SBarry Smith   PetscInt       i;
14698c34d3f5SBarry Smith   PetscOption   amsopt;
147053acd3b1SBarry Smith 
147153acd3b1SBarry Smith   PetscFunctionBegin;
1472e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1473e26ddf31SBarry Smith     PetscInt *vals;
1474e26ddf31SBarry Smith 
1475e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1476854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1477e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1478e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1479e26ddf31SBarry Smith     amsopt->arraylength = *n;
1480e26ddf31SBarry Smith   }
1481e55864a3SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1482e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1483e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
148453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1485e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
148653acd3b1SBarry Smith     }
1487e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
148853acd3b1SBarry Smith   }
148953acd3b1SBarry Smith   PetscFunctionReturn(0);
149053acd3b1SBarry Smith }
149153acd3b1SBarry Smith 
149253acd3b1SBarry Smith #undef __FUNCT__
1493e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
149453acd3b1SBarry Smith /*@C
149553acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
149653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
149753acd3b1SBarry Smith    no intervening spaces.
149853acd3b1SBarry Smith 
14993f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
150053acd3b1SBarry Smith 
150153acd3b1SBarry Smith    Input Parameters:
150253acd3b1SBarry Smith +  opt - the option one is seeking
150353acd3b1SBarry Smith .  text - short string describing option
150453acd3b1SBarry Smith .  man - manual page for option
150553acd3b1SBarry Smith -  nmax - maximum number of strings
150653acd3b1SBarry Smith 
150753acd3b1SBarry Smith    Output Parameter:
150853acd3b1SBarry Smith +  value - location to copy strings
150953acd3b1SBarry Smith .  nmax - actual number of strings found
151053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
151153acd3b1SBarry Smith 
151253acd3b1SBarry Smith    Level: beginner
151353acd3b1SBarry Smith 
151453acd3b1SBarry Smith    Notes:
151553acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
151653acd3b1SBarry Smith    strings returned by this function.
151753acd3b1SBarry Smith 
151853acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
151953acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
152053acd3b1SBarry Smith 
152153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
152253acd3b1SBarry Smith 
152353acd3b1SBarry Smith    Concepts: options database^array of strings
152453acd3b1SBarry Smith 
152553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1526acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
152753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
152853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1529acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1530a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
153153acd3b1SBarry Smith @*/
15328c34d3f5SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
153353acd3b1SBarry Smith {
153453acd3b1SBarry Smith   PetscErrorCode ierr;
15358c34d3f5SBarry Smith   PetscOption   amsopt;
153653acd3b1SBarry Smith 
153753acd3b1SBarry Smith   PetscFunctionBegin;
1538e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1539e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1540854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1541a297a907SKarl Rupp 
15421ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
15431ae3d29cSBarry Smith   }
1544e55864a3SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1545e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1546e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
154753acd3b1SBarry Smith   }
154853acd3b1SBarry Smith   PetscFunctionReturn(0);
154953acd3b1SBarry Smith }
155053acd3b1SBarry Smith 
1551e2446a98SMatthew Knepley #undef __FUNCT__
1552e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1553e2446a98SMatthew Knepley /*@C
1554acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1555e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1556e2446a98SMatthew Knepley    no intervening spaces.
1557e2446a98SMatthew Knepley 
15583f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1559e2446a98SMatthew Knepley 
1560e2446a98SMatthew Knepley    Input Parameters:
1561e2446a98SMatthew Knepley +  opt - the option one is seeking
1562e2446a98SMatthew Knepley .  text - short string describing option
1563e2446a98SMatthew Knepley .  man - manual page for option
1564e2446a98SMatthew Knepley -  nmax - maximum number of values
1565e2446a98SMatthew Knepley 
1566e2446a98SMatthew Knepley    Output Parameter:
1567e2446a98SMatthew Knepley +  value - location to copy values
1568e2446a98SMatthew Knepley .  nmax - actual number of values found
1569e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1570e2446a98SMatthew Knepley 
1571e2446a98SMatthew Knepley    Level: beginner
1572e2446a98SMatthew Knepley 
1573e2446a98SMatthew Knepley    Notes:
1574e2446a98SMatthew Knepley    The user should pass in an array of doubles
1575e2446a98SMatthew Knepley 
1576e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1577e2446a98SMatthew Knepley 
1578e2446a98SMatthew Knepley    Concepts: options database^array of strings
1579e2446a98SMatthew Knepley 
1580e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1581acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1582e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1583e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1584acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1585a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1586e2446a98SMatthew Knepley @*/
15878c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1588e2446a98SMatthew Knepley {
1589e2446a98SMatthew Knepley   PetscErrorCode ierr;
1590e2446a98SMatthew Knepley   PetscInt       i;
15918c34d3f5SBarry Smith   PetscOption    amsopt;
1592e2446a98SMatthew Knepley 
1593e2446a98SMatthew Knepley   PetscFunctionBegin;
1594e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1595ace3abfcSBarry Smith     PetscBool *vals;
15961ae3d29cSBarry Smith 
159783355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
15981a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1599ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
16001ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
16011ae3d29cSBarry Smith     amsopt->arraylength = *n;
16021ae3d29cSBarry Smith   }
1603e55864a3SBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1604e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1605e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1606e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1607e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1608e2446a98SMatthew Knepley     }
1609e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1610e2446a98SMatthew Knepley   }
1611e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1612e2446a98SMatthew Knepley }
1613e2446a98SMatthew Knepley 
16148cc676e6SMatthew G Knepley #undef __FUNCT__
1615e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
16168cc676e6SMatthew G Knepley /*@C
1617d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
16188cc676e6SMatthew G Knepley 
16198cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
16208cc676e6SMatthew G Knepley 
16218cc676e6SMatthew G Knepley    Input Parameters:
16228cc676e6SMatthew G Knepley +  opt - option name
16238cc676e6SMatthew G Knepley .  text - short string that describes the option
16248cc676e6SMatthew G Knepley -  man - manual page with additional information on option
16258cc676e6SMatthew G Knepley 
16268cc676e6SMatthew G Knepley    Output Parameter:
16278cc676e6SMatthew G Knepley +  viewer - the viewer
16288cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
16298cc676e6SMatthew G Knepley 
16308cc676e6SMatthew G Knepley    Level: beginner
16318cc676e6SMatthew G Knepley 
16328cc676e6SMatthew G Knepley    Concepts: options database^has int
16338cc676e6SMatthew G Knepley 
16348cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
16358cc676e6SMatthew G Knepley 
16365a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
16378cc676e6SMatthew G Knepley 
16388cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
16398cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
16408cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
16418cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
16428cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
16438cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1644a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
16458cc676e6SMatthew G Knepley @*/
16468c34d3f5SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
16478cc676e6SMatthew G Knepley {
16488cc676e6SMatthew G Knepley   PetscErrorCode ierr;
16498c34d3f5SBarry Smith   PetscOption    amsopt;
16508cc676e6SMatthew G Knepley 
16518cc676e6SMatthew G Knepley   PetscFunctionBegin;
16521a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
16531a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
165464facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
16555b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
16568cc676e6SMatthew G Knepley   }
1657e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1658e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1659e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
16608cc676e6SMatthew G Knepley   }
16618cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
16628cc676e6SMatthew G Knepley }
16638cc676e6SMatthew G Knepley 
166453acd3b1SBarry Smith 
166553acd3b1SBarry Smith #undef __FUNCT__
166653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
166753acd3b1SBarry Smith /*@C
1668b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
166953acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
167053acd3b1SBarry Smith 
16713f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
167253acd3b1SBarry Smith 
167353acd3b1SBarry Smith    Input Parameter:
167453acd3b1SBarry Smith .   head - the heading text
167553acd3b1SBarry Smith 
167653acd3b1SBarry Smith 
167753acd3b1SBarry Smith    Level: intermediate
167853acd3b1SBarry Smith 
167953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
168053acd3b1SBarry Smith 
1681b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
168253acd3b1SBarry Smith 
168353acd3b1SBarry Smith    Concepts: options database^subheading
168453acd3b1SBarry Smith 
168553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1686acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
168753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
168853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1689acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1690a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
169153acd3b1SBarry Smith @*/
16928c34d3f5SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptions *PetscOptionsObject,const char head[])
169353acd3b1SBarry Smith {
169453acd3b1SBarry Smith   PetscErrorCode ierr;
169553acd3b1SBarry Smith 
169653acd3b1SBarry Smith   PetscFunctionBegin;
1697e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1698e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
169953acd3b1SBarry Smith   }
170053acd3b1SBarry Smith   PetscFunctionReturn(0);
170153acd3b1SBarry Smith }
170253acd3b1SBarry Smith 
170353acd3b1SBarry Smith 
170453acd3b1SBarry Smith 
170553acd3b1SBarry Smith 
170653acd3b1SBarry Smith 
170753acd3b1SBarry Smith 
1708