xref: /petsc/src/sys/objects/aoptions.c (revision 6991f827c8a892fc5b642fee1563a9b520cca98c)
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;
5857781c08eSBarry Smith       case OPTION_BOOL:
586e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5876356e834SBarry Smith         break;
5887781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
589e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
590e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
591e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5921ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5931ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5941ae3d29cSBarry Smith         }
5951ae3d29cSBarry Smith         break;
596a264d7a6SBarry Smith       case OPTION_FLIST:
597*6991f827SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
598*6991f827SBarry Smith         break;
5991ae3d29cSBarry Smith       case OPTION_ELIST:
600e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
6016356e834SBarry Smith         break;
6021ae3d29cSBarry Smith       case OPTION_STRING:
603e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
604d51da6bfSBarry Smith         break;
6051ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
606e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
607e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
608e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
6091ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
6101ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
6111ae3d29cSBarry Smith         }
6126356e834SBarry Smith         break;
6136356e834SBarry Smith       }
6146356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
6156356e834SBarry Smith     }
616*6991f827SBarry Smith     if (PetscOptionsObject->next->type == OPTION_ELIST) {
617*6991f827SBarry Smith       ierr = PetscStrNArrayDestroy(PetscOptionsObject->next->nlist,(char ***)&PetscOptionsObject->next->list);CHKERRQ(ierr);
618*6991f827SBarry Smith     }
619e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
620e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
621e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
622e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
623c979a496SBarry Smith 
62483355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
62583355fc5SBarry Smith       free(PetscOptionsObject->next->data);
626c979a496SBarry Smith     } else {
62783355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
628c979a496SBarry Smith     }
6297781c08eSBarry Smith 
63083355fc5SBarry Smith     last                    = PetscOptionsObject->next;
63183355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6326356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6336356e834SBarry Smith   }
634f59f755dSBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
635e55864a3SBarry Smith   PetscOptionsObject->next = 0;
63653acd3b1SBarry Smith   PetscFunctionReturn(0);
63753acd3b1SBarry Smith }
63853acd3b1SBarry Smith 
63953acd3b1SBarry Smith #undef __FUNCT__
640e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
64153acd3b1SBarry Smith /*@C
64253acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
64353acd3b1SBarry Smith 
6443f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
64553acd3b1SBarry Smith 
64653acd3b1SBarry Smith    Input Parameters:
64753acd3b1SBarry Smith +  opt - option name
64853acd3b1SBarry Smith .  text - short string that describes the option
64953acd3b1SBarry Smith .  man - manual page with additional information on option
65053acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6510fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6520fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6530fdccdaeSBarry Smith $                 value = defaultvalue
6540fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6550fdccdaeSBarry Smith $                 if (flg) {
65653acd3b1SBarry Smith 
65753acd3b1SBarry Smith    Output Parameter:
65853acd3b1SBarry Smith +  value - the  value to return
659b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
66053acd3b1SBarry Smith 
66153acd3b1SBarry Smith    Level: beginner
66253acd3b1SBarry Smith 
66353acd3b1SBarry Smith    Concepts: options database
66453acd3b1SBarry Smith 
66553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
66653acd3b1SBarry Smith 
66753acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
66853acd3b1SBarry Smith 
66953acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
670acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
671acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
67253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
67353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
674acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
675a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
67653acd3b1SBarry Smith @*/
6778c34d3f5SBarry 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)
67853acd3b1SBarry Smith {
67953acd3b1SBarry Smith   PetscErrorCode ierr;
68053acd3b1SBarry Smith   PetscInt       ntext = 0;
681aa5bb8c0SSatish Balay   PetscInt       tval;
682ace3abfcSBarry Smith   PetscBool      tflg;
68353acd3b1SBarry Smith 
68453acd3b1SBarry Smith   PetscFunctionBegin;
68553acd3b1SBarry Smith   while (list[ntext++]) {
686e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
68753acd3b1SBarry Smith   }
688e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
68953acd3b1SBarry Smith   ntext -= 3;
690e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
691aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
692aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
693aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
69453acd3b1SBarry Smith   PetscFunctionReturn(0);
69553acd3b1SBarry Smith }
69653acd3b1SBarry Smith 
69753acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
69853acd3b1SBarry Smith #undef __FUNCT__
699e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
70053acd3b1SBarry Smith /*@C
70153acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
70253acd3b1SBarry Smith 
7033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith    Input Parameters:
70653acd3b1SBarry Smith +  opt - option name
70753acd3b1SBarry Smith .  text - short string that describes the option
70853acd3b1SBarry Smith .  man - manual page with additional information on option
7090fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7100fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
7110fdccdaeSBarry Smith $                 value = defaultvalue
7120fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
7130fdccdaeSBarry Smith $                 if (flg) {
71453acd3b1SBarry Smith 
71553acd3b1SBarry Smith    Output Parameter:
71653acd3b1SBarry Smith +  value - the integer value to return
71753acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
71853acd3b1SBarry Smith 
71953acd3b1SBarry Smith    Level: beginner
72053acd3b1SBarry Smith 
72153acd3b1SBarry Smith    Concepts: options database^has int
72253acd3b1SBarry Smith 
72353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
72453acd3b1SBarry Smith 
72553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
726acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
727acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
72853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
72953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
730acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
731a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
73253acd3b1SBarry Smith @*/
7338c34d3f5SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
73453acd3b1SBarry Smith {
73553acd3b1SBarry Smith   PetscErrorCode ierr;
7368c34d3f5SBarry Smith   PetscOption    amsopt;
73712655325SBarry Smith   PetscBool      wasset;
73853acd3b1SBarry Smith 
73953acd3b1SBarry Smith   PetscFunctionBegin;
740e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
741e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7426356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
74312655325SBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
7443e211508SBarry Smith 
74512655325SBarry Smith     ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,&currentvalue,&wasset);CHKERRQ(ierr);
7463e211508SBarry Smith     if (wasset) {
74712655325SBarry Smith       *(PetscInt*)amsopt->data = currentvalue;
7483e211508SBarry Smith     }
749af6d86caSBarry Smith   }
750e55864a3SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
751e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
7521a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
75353acd3b1SBarry Smith   }
75453acd3b1SBarry Smith   PetscFunctionReturn(0);
75553acd3b1SBarry Smith }
75653acd3b1SBarry Smith 
75753acd3b1SBarry Smith #undef __FUNCT__
758e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
75953acd3b1SBarry Smith /*@C
76053acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
76153acd3b1SBarry Smith 
7623f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
76353acd3b1SBarry Smith 
76453acd3b1SBarry Smith    Input Parameters:
76553acd3b1SBarry Smith +  opt - option name
76653acd3b1SBarry Smith .  text - short string that describes the option
76753acd3b1SBarry Smith .  man - manual page with additional information on option
7680fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
769bcbf2dc5SJed Brown -  len - length of the result string including null terminator
77053acd3b1SBarry Smith 
77153acd3b1SBarry Smith    Output Parameter:
77253acd3b1SBarry Smith +  value - the value to return
77353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
77453acd3b1SBarry Smith 
77553acd3b1SBarry Smith    Level: beginner
77653acd3b1SBarry Smith 
77753acd3b1SBarry Smith    Concepts: options database^has int
77853acd3b1SBarry Smith 
77953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
78053acd3b1SBarry Smith 
7817fccdfe4SBarry 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).
7827fccdfe4SBarry Smith 
78353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
784acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
785acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
78653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
78753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
788acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
789a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
79053acd3b1SBarry Smith @*/
7918c34d3f5SBarry 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)
79253acd3b1SBarry Smith {
79353acd3b1SBarry Smith   PetscErrorCode ierr;
7948c34d3f5SBarry Smith   PetscOption   amsopt;
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith   PetscFunctionBegin;
7971a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
7981a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
79964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
8000fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
801af6d86caSBarry Smith   }
802e55864a3SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
803e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8041a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
80553acd3b1SBarry Smith   }
80653acd3b1SBarry Smith   PetscFunctionReturn(0);
80753acd3b1SBarry Smith }
80853acd3b1SBarry Smith 
80953acd3b1SBarry Smith #undef __FUNCT__
810e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
81153acd3b1SBarry Smith /*@C
81253acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
81353acd3b1SBarry Smith 
8143f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
81553acd3b1SBarry Smith 
81653acd3b1SBarry Smith    Input Parameters:
81753acd3b1SBarry Smith +  opt - option name
81853acd3b1SBarry Smith .  text - short string that describes the option
81953acd3b1SBarry Smith .  man - manual page with additional information on option
8200fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8210fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
8220fdccdaeSBarry Smith $                 value = defaultvalue
8230fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
8240fdccdaeSBarry Smith $                 if (flg) {
82553acd3b1SBarry Smith 
82653acd3b1SBarry Smith    Output Parameter:
82753acd3b1SBarry Smith +  value - the value to return
82853acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
82953acd3b1SBarry Smith 
83053acd3b1SBarry Smith    Level: beginner
83153acd3b1SBarry Smith 
83253acd3b1SBarry Smith    Concepts: options database^has int
83353acd3b1SBarry Smith 
83453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
83553acd3b1SBarry Smith 
83653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
837acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
838acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
83953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
84053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
841acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
842a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
84353acd3b1SBarry Smith @*/
8448c34d3f5SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
84553acd3b1SBarry Smith {
84653acd3b1SBarry Smith   PetscErrorCode ierr;
8478c34d3f5SBarry Smith   PetscOption   amsopt;
84853acd3b1SBarry Smith 
84953acd3b1SBarry Smith   PetscFunctionBegin;
850e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
851e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
852538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
853a297a907SKarl Rupp 
8540fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
855538aa990SBarry Smith   }
8561a1499c8SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
8571a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
8581a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
85953acd3b1SBarry Smith   }
86053acd3b1SBarry Smith   PetscFunctionReturn(0);
86153acd3b1SBarry Smith }
86253acd3b1SBarry Smith 
86353acd3b1SBarry Smith #undef __FUNCT__
864e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
86553acd3b1SBarry Smith /*@C
86653acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
86753acd3b1SBarry Smith 
8683f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
86953acd3b1SBarry Smith 
87053acd3b1SBarry Smith    Input Parameters:
87153acd3b1SBarry Smith +  opt - option name
87253acd3b1SBarry Smith .  text - short string that describes the option
87353acd3b1SBarry Smith .  man - manual page with additional information on option
8740fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8750fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
8760fdccdaeSBarry Smith $                 value = defaultvalue
8770fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
8780fdccdaeSBarry Smith $                 if (flg) {
8790fdccdaeSBarry Smith 
88053acd3b1SBarry Smith 
88153acd3b1SBarry Smith    Output Parameter:
88253acd3b1SBarry Smith +  value - the value to return
88353acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
88453acd3b1SBarry Smith 
88553acd3b1SBarry Smith    Level: beginner
88653acd3b1SBarry Smith 
88753acd3b1SBarry Smith    Concepts: options database^has int
88853acd3b1SBarry Smith 
88953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89053acd3b1SBarry Smith 
89153acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
892acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
893acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
89453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
89553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
896acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
897a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
89853acd3b1SBarry Smith @*/
8998c34d3f5SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
90053acd3b1SBarry Smith {
90153acd3b1SBarry Smith   PetscErrorCode ierr;
90253acd3b1SBarry Smith 
90353acd3b1SBarry Smith   PetscFunctionBegin;
90453acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
9050fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
90653acd3b1SBarry Smith #else
907e55864a3SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
90853acd3b1SBarry Smith #endif
90953acd3b1SBarry Smith   PetscFunctionReturn(0);
91053acd3b1SBarry Smith }
91153acd3b1SBarry Smith 
91253acd3b1SBarry Smith #undef __FUNCT__
913e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
91453acd3b1SBarry Smith /*@C
91590d69ab7SBarry 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
91690d69ab7SBarry Smith                       its value is set to false.
91753acd3b1SBarry Smith 
9183f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
91953acd3b1SBarry Smith 
92053acd3b1SBarry Smith    Input Parameters:
92153acd3b1SBarry Smith +  opt - option name
92253acd3b1SBarry Smith .  text - short string that describes the option
92353acd3b1SBarry Smith -  man - manual page with additional information on option
92453acd3b1SBarry Smith 
92553acd3b1SBarry Smith    Output Parameter:
92653acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
92753acd3b1SBarry Smith 
92853acd3b1SBarry Smith    Level: beginner
92953acd3b1SBarry Smith 
93053acd3b1SBarry Smith    Concepts: options database^has int
93153acd3b1SBarry Smith 
93253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
93353acd3b1SBarry Smith 
93453acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
935acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
936acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
93753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
93853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
939acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
940a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
94153acd3b1SBarry Smith @*/
9428c34d3f5SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
94353acd3b1SBarry Smith {
94453acd3b1SBarry Smith   PetscErrorCode ierr;
9458c34d3f5SBarry Smith   PetscOption   amsopt;
94653acd3b1SBarry Smith 
94753acd3b1SBarry Smith   PetscFunctionBegin;
948e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
949e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
950ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
951a297a907SKarl Rupp 
952ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9531ae3d29cSBarry Smith   }
954e55864a3SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
955e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
956e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
95753acd3b1SBarry Smith   }
95853acd3b1SBarry Smith   PetscFunctionReturn(0);
95953acd3b1SBarry Smith }
96053acd3b1SBarry Smith 
96153acd3b1SBarry Smith #undef __FUNCT__
962e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
96353acd3b1SBarry Smith /*@C
964a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
96553acd3b1SBarry Smith 
9663f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
96753acd3b1SBarry Smith 
96853acd3b1SBarry Smith    Input Parameters:
96953acd3b1SBarry Smith +  opt - option name
97053acd3b1SBarry Smith .  text - short string that describes the option
97153acd3b1SBarry Smith .  man - manual page with additional information on option
97253acd3b1SBarry Smith .  list - the possible choices
9730fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
9740fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
9750fdccdaeSBarry Smith $                 if (flg) {
9763cc1e11dSBarry Smith -  len - the length of the character array value
97753acd3b1SBarry Smith 
97853acd3b1SBarry Smith    Output Parameter:
97953acd3b1SBarry Smith +  value - the value to return
98053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
98153acd3b1SBarry Smith 
98253acd3b1SBarry Smith    Level: intermediate
98353acd3b1SBarry Smith 
98453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
98553acd3b1SBarry Smith 
98653acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
98753acd3b1SBarry Smith 
98853acd3b1SBarry Smith    To get a listing of all currently specified options,
98988c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
99053acd3b1SBarry Smith 
991eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
992eabe10d7SBarry Smith 
99353acd3b1SBarry Smith    Concepts: options database^list
99453acd3b1SBarry Smith 
99553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
996acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
99753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
99853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
999acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1000a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
100153acd3b1SBarry Smith @*/
10028c34d3f5SBarry 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)
100353acd3b1SBarry Smith {
100453acd3b1SBarry Smith   PetscErrorCode ierr;
10058c34d3f5SBarry Smith   PetscOption   amsopt;
100653acd3b1SBarry Smith 
100753acd3b1SBarry Smith   PetscFunctionBegin;
10081a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10091a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
101064facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10110fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10123cc1e11dSBarry Smith     amsopt->flist = list;
10133cc1e11dSBarry Smith   }
10141a1499c8SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
10151a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10161a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
101753acd3b1SBarry Smith   }
101853acd3b1SBarry Smith   PetscFunctionReturn(0);
101953acd3b1SBarry Smith }
102053acd3b1SBarry Smith 
102153acd3b1SBarry Smith #undef __FUNCT__
1022e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
102353acd3b1SBarry Smith /*@C
102453acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
102553acd3b1SBarry Smith 
10263f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
102753acd3b1SBarry Smith 
102853acd3b1SBarry Smith    Input Parameters:
102953acd3b1SBarry Smith +  opt - option name
103053acd3b1SBarry Smith .  ltext - short string that describes the option
103153acd3b1SBarry Smith .  man - manual page with additional information on option
1032a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
103353acd3b1SBarry Smith .  ntext - number of choices
10340fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
10350fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
10360fdccdaeSBarry Smith $                 if (flg) {
10370fdccdaeSBarry Smith 
103853acd3b1SBarry Smith 
103953acd3b1SBarry Smith    Output Parameter:
104053acd3b1SBarry Smith +  value - the index of the value to return
104153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
104253acd3b1SBarry Smith 
104353acd3b1SBarry Smith    Level: intermediate
104453acd3b1SBarry Smith 
104553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
104653acd3b1SBarry Smith 
1047a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
104853acd3b1SBarry Smith 
104953acd3b1SBarry Smith    Concepts: options database^list
105053acd3b1SBarry Smith 
105153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1052acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
105353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
105453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1055acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1056a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
105753acd3b1SBarry Smith @*/
10588c34d3f5SBarry 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)
105953acd3b1SBarry Smith {
106053acd3b1SBarry Smith   PetscErrorCode ierr;
106153acd3b1SBarry Smith   PetscInt       i;
10628c34d3f5SBarry Smith   PetscOption   amsopt;
106353acd3b1SBarry Smith 
106453acd3b1SBarry Smith   PetscFunctionBegin;
10651a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
10661a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
106764facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10680fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
1069*6991f827SBarry Smith     ierr = PetscStrNArrayallocpy(ntext,list,(char***)&amsopt->list);CHKERRQ(ierr);
10701ae3d29cSBarry Smith     amsopt->nlist = ntext;
10711ae3d29cSBarry Smith   }
10721a1499c8SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
10731a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
10741a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
107553acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1076e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
107753acd3b1SBarry Smith     }
1078e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
107953acd3b1SBarry Smith   }
108053acd3b1SBarry Smith   PetscFunctionReturn(0);
108153acd3b1SBarry Smith }
108253acd3b1SBarry Smith 
108353acd3b1SBarry Smith #undef __FUNCT__
1084e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
108553acd3b1SBarry Smith /*@C
1086acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1087d5649816SBarry Smith        which at most a single value can be true.
108853acd3b1SBarry Smith 
10893f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
109053acd3b1SBarry Smith 
109153acd3b1SBarry Smith    Input Parameters:
109253acd3b1SBarry Smith +  opt - option name
109353acd3b1SBarry Smith .  text - short string that describes the option
109453acd3b1SBarry Smith -  man - manual page with additional information on option
109553acd3b1SBarry Smith 
109653acd3b1SBarry Smith    Output Parameter:
109753acd3b1SBarry Smith .  flg - whether that option was set or not
109853acd3b1SBarry Smith 
109953acd3b1SBarry Smith    Level: intermediate
110053acd3b1SBarry Smith 
110153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
110253acd3b1SBarry Smith 
1103acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
110453acd3b1SBarry Smith 
110553acd3b1SBarry Smith     Concepts: options database^logical group
110653acd3b1SBarry Smith 
110753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1108acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
110953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
111053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1111acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1112a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
111353acd3b1SBarry Smith @*/
11148c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
111553acd3b1SBarry Smith {
111653acd3b1SBarry Smith   PetscErrorCode ierr;
11178c34d3f5SBarry Smith   PetscOption   amsopt;
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith   PetscFunctionBegin;
1120e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
112183355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1122ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1123a297a907SKarl Rupp 
1124ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11251ae3d29cSBarry Smith   }
112668b16fdaSBarry Smith   *flg = PETSC_FALSE;
1127e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1128e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1129e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1130e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
113153acd3b1SBarry Smith   }
113253acd3b1SBarry Smith   PetscFunctionReturn(0);
113353acd3b1SBarry Smith }
113453acd3b1SBarry Smith 
113553acd3b1SBarry Smith #undef __FUNCT__
1136e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
113753acd3b1SBarry Smith /*@C
1138acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1139d5649816SBarry Smith        which at most a single value can be true.
114053acd3b1SBarry Smith 
11413f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
114253acd3b1SBarry Smith 
114353acd3b1SBarry Smith    Input Parameters:
114453acd3b1SBarry Smith +  opt - option name
114553acd3b1SBarry Smith .  text - short string that describes the option
114653acd3b1SBarry Smith -  man - manual page with additional information on option
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith    Output Parameter:
114953acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
115053acd3b1SBarry Smith 
115153acd3b1SBarry Smith    Level: intermediate
115253acd3b1SBarry Smith 
115353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
115453acd3b1SBarry Smith 
1155acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
115653acd3b1SBarry Smith 
115753acd3b1SBarry Smith     Concepts: options database^logical group
115853acd3b1SBarry Smith 
115953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1160acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
116153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
116253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1163acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1164a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
116553acd3b1SBarry Smith @*/
11668c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
116753acd3b1SBarry Smith {
116853acd3b1SBarry Smith   PetscErrorCode ierr;
11698c34d3f5SBarry Smith   PetscOption   amsopt;
117053acd3b1SBarry Smith 
117153acd3b1SBarry Smith   PetscFunctionBegin;
1172e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
117383355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1174ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1175a297a907SKarl Rupp 
1176ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11771ae3d29cSBarry Smith   }
117817326d04SJed Brown   *flg = PETSC_FALSE;
1179e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1180e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1181e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
118253acd3b1SBarry Smith   }
118353acd3b1SBarry Smith   PetscFunctionReturn(0);
118453acd3b1SBarry Smith }
118553acd3b1SBarry Smith 
118653acd3b1SBarry Smith #undef __FUNCT__
1187e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
118853acd3b1SBarry Smith /*@C
1189acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1190d5649816SBarry Smith        which at most a single value can be true.
119153acd3b1SBarry Smith 
11923f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
119353acd3b1SBarry Smith 
119453acd3b1SBarry Smith    Input Parameters:
119553acd3b1SBarry Smith +  opt - option name
119653acd3b1SBarry Smith .  text - short string that describes the option
119753acd3b1SBarry Smith -  man - manual page with additional information on option
119853acd3b1SBarry Smith 
119953acd3b1SBarry Smith    Output Parameter:
120053acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
120153acd3b1SBarry Smith 
120253acd3b1SBarry Smith    Level: intermediate
120353acd3b1SBarry Smith 
120453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
120553acd3b1SBarry Smith 
1206acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
120753acd3b1SBarry Smith 
120853acd3b1SBarry Smith     Concepts: options database^logical group
120953acd3b1SBarry Smith 
121053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1211acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
121253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
121353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1214acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1215a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
121653acd3b1SBarry Smith @*/
12178c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
121853acd3b1SBarry Smith {
121953acd3b1SBarry Smith   PetscErrorCode ierr;
12208c34d3f5SBarry Smith   PetscOption   amsopt;
122153acd3b1SBarry Smith 
122253acd3b1SBarry Smith   PetscFunctionBegin;
1223e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
122483355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1225ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1226a297a907SKarl Rupp 
1227ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
12281ae3d29cSBarry Smith   }
122917326d04SJed Brown   *flg = PETSC_FALSE;
1230e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1231e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1232e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
123353acd3b1SBarry Smith   }
123453acd3b1SBarry Smith   PetscFunctionReturn(0);
123553acd3b1SBarry Smith }
123653acd3b1SBarry Smith 
123753acd3b1SBarry Smith #undef __FUNCT__
1238e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
123953acd3b1SBarry Smith /*@C
1240acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
124153acd3b1SBarry Smith 
12423f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
124353acd3b1SBarry Smith 
124453acd3b1SBarry Smith    Input Parameters:
124553acd3b1SBarry Smith +  opt - option name
124653acd3b1SBarry Smith .  text - short string that describes the option
1247868c398cSBarry Smith .  man - manual page with additional information on option
124894ae4db5SBarry Smith -  currentvalue - the current value
124953acd3b1SBarry Smith 
125053acd3b1SBarry Smith    Output Parameter:
125153acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
125253acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
125353acd3b1SBarry Smith 
125453acd3b1SBarry Smith    Level: beginner
125553acd3b1SBarry Smith 
125653acd3b1SBarry Smith    Concepts: options database^logical
125753acd3b1SBarry Smith 
125853acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
125953acd3b1SBarry Smith 
126053acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1261acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1262acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
126353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
126453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1265acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1266a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
126753acd3b1SBarry Smith @*/
12688c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
126953acd3b1SBarry Smith {
127053acd3b1SBarry Smith   PetscErrorCode ierr;
1271ace3abfcSBarry Smith   PetscBool      iset;
12728c34d3f5SBarry Smith   PetscOption   amsopt;
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith   PetscFunctionBegin;
1275e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1276e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1277ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1278a297a907SKarl Rupp 
127994ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1280af6d86caSBarry Smith   }
12811a1499c8SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
128253acd3b1SBarry Smith   if (set) *set = iset;
12831a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
128494ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
12851a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
128653acd3b1SBarry Smith   }
128753acd3b1SBarry Smith   PetscFunctionReturn(0);
128853acd3b1SBarry Smith }
128953acd3b1SBarry Smith 
129053acd3b1SBarry Smith #undef __FUNCT__
1291e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
129253acd3b1SBarry Smith /*@C
129353acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
129453acd3b1SBarry Smith    option in the database. The values must be separated with commas with
129553acd3b1SBarry Smith    no intervening spaces.
129653acd3b1SBarry Smith 
12973f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
129853acd3b1SBarry Smith 
129953acd3b1SBarry Smith    Input Parameters:
130053acd3b1SBarry Smith +  opt - the option one is seeking
130153acd3b1SBarry Smith .  text - short string describing option
130253acd3b1SBarry Smith .  man - manual page for option
130353acd3b1SBarry Smith -  nmax - maximum number of values
130453acd3b1SBarry Smith 
130553acd3b1SBarry Smith    Output Parameter:
130653acd3b1SBarry Smith +  value - location to copy values
130753acd3b1SBarry Smith .  nmax - actual number of values found
130853acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
130953acd3b1SBarry Smith 
131053acd3b1SBarry Smith    Level: beginner
131153acd3b1SBarry Smith 
131253acd3b1SBarry Smith    Notes:
131353acd3b1SBarry Smith    The user should pass in an array of doubles
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
131653acd3b1SBarry Smith 
131753acd3b1SBarry Smith    Concepts: options database^array of strings
131853acd3b1SBarry Smith 
131953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1320acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
132153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
132253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1323acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1324a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
132553acd3b1SBarry Smith @*/
13268c34d3f5SBarry Smith PetscErrorCode  PetscOptionsRealArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
132753acd3b1SBarry Smith {
132853acd3b1SBarry Smith   PetscErrorCode ierr;
132953acd3b1SBarry Smith   PetscInt       i;
13308c34d3f5SBarry Smith   PetscOption   amsopt;
133153acd3b1SBarry Smith 
133253acd3b1SBarry Smith   PetscFunctionBegin;
1333e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1334e26ddf31SBarry Smith     PetscReal *vals;
1335e26ddf31SBarry Smith 
1336e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1337e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1338e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1339e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1340e26ddf31SBarry Smith     amsopt->arraylength = *n;
1341e26ddf31SBarry Smith   }
1342e55864a3SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1343e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1344a519f713SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,(double)value[0]);CHKERRQ(ierr);
134553acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1346e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
134753acd3b1SBarry Smith     }
1348e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
134953acd3b1SBarry Smith   }
135053acd3b1SBarry Smith   PetscFunctionReturn(0);
135153acd3b1SBarry Smith }
135253acd3b1SBarry Smith 
135353acd3b1SBarry Smith 
135453acd3b1SBarry Smith #undef __FUNCT__
1355e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
135653acd3b1SBarry Smith /*@C
135753acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1358b32a342fSShri Abhyankar    option in the database.
135953acd3b1SBarry Smith 
13603f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
136153acd3b1SBarry Smith 
136253acd3b1SBarry Smith    Input Parameters:
136353acd3b1SBarry Smith +  opt - the option one is seeking
136453acd3b1SBarry Smith .  text - short string describing option
136553acd3b1SBarry Smith .  man - manual page for option
1366f8a50e2bSBarry Smith -  n - maximum number of values
136753acd3b1SBarry Smith 
136853acd3b1SBarry Smith    Output Parameter:
136953acd3b1SBarry Smith +  value - location to copy values
1370f8a50e2bSBarry Smith .  n - actual number of values found
137153acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
137253acd3b1SBarry Smith 
137353acd3b1SBarry Smith    Level: beginner
137453acd3b1SBarry Smith 
137553acd3b1SBarry Smith    Notes:
1376b32a342fSShri Abhyankar    The array can be passed as
1377b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13780fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13790fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13800fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1381b32a342fSShri Abhyankar 
1382b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
138353acd3b1SBarry Smith 
138453acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
138553acd3b1SBarry Smith 
1386b32a342fSShri Abhyankar    Concepts: options database^array of ints
138753acd3b1SBarry Smith 
138853acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1389acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
139053acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
139153acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1392acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1393a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
139453acd3b1SBarry Smith @*/
13958c34d3f5SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
139653acd3b1SBarry Smith {
139753acd3b1SBarry Smith   PetscErrorCode ierr;
139853acd3b1SBarry Smith   PetscInt       i;
13998c34d3f5SBarry Smith   PetscOption   amsopt;
140053acd3b1SBarry Smith 
140153acd3b1SBarry Smith   PetscFunctionBegin;
1402e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1403e26ddf31SBarry Smith     PetscInt *vals;
1404e26ddf31SBarry Smith 
1405e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1406854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1407e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1408e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1409e26ddf31SBarry Smith     amsopt->arraylength = *n;
1410e26ddf31SBarry Smith   }
1411e55864a3SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1412e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1413e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
141453acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1415e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
141653acd3b1SBarry Smith     }
1417e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
141853acd3b1SBarry Smith   }
141953acd3b1SBarry Smith   PetscFunctionReturn(0);
142053acd3b1SBarry Smith }
142153acd3b1SBarry Smith 
142253acd3b1SBarry Smith #undef __FUNCT__
1423e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
142453acd3b1SBarry Smith /*@C
142553acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
142653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
142753acd3b1SBarry Smith    no intervening spaces.
142853acd3b1SBarry Smith 
14293f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
143053acd3b1SBarry Smith 
143153acd3b1SBarry Smith    Input Parameters:
143253acd3b1SBarry Smith +  opt - the option one is seeking
143353acd3b1SBarry Smith .  text - short string describing option
143453acd3b1SBarry Smith .  man - manual page for option
143553acd3b1SBarry Smith -  nmax - maximum number of strings
143653acd3b1SBarry Smith 
143753acd3b1SBarry Smith    Output Parameter:
143853acd3b1SBarry Smith +  value - location to copy strings
143953acd3b1SBarry Smith .  nmax - actual number of strings found
144053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
144153acd3b1SBarry Smith 
144253acd3b1SBarry Smith    Level: beginner
144353acd3b1SBarry Smith 
144453acd3b1SBarry Smith    Notes:
144553acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
144653acd3b1SBarry Smith    strings returned by this function.
144753acd3b1SBarry Smith 
144853acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
144953acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
145053acd3b1SBarry Smith 
145153acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
145253acd3b1SBarry Smith 
145353acd3b1SBarry Smith    Concepts: options database^array of strings
145453acd3b1SBarry Smith 
145553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1456acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
145753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
145853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1459acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1460a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
146153acd3b1SBarry Smith @*/
14628c34d3f5SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
146353acd3b1SBarry Smith {
146453acd3b1SBarry Smith   PetscErrorCode ierr;
14658c34d3f5SBarry Smith   PetscOption   amsopt;
146653acd3b1SBarry Smith 
146753acd3b1SBarry Smith   PetscFunctionBegin;
1468e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1469e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1470854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1471a297a907SKarl Rupp 
14721ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14731ae3d29cSBarry Smith   }
1474e55864a3SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1475e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1476e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
147753acd3b1SBarry Smith   }
147853acd3b1SBarry Smith   PetscFunctionReturn(0);
147953acd3b1SBarry Smith }
148053acd3b1SBarry Smith 
1481e2446a98SMatthew Knepley #undef __FUNCT__
1482e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1483e2446a98SMatthew Knepley /*@C
1484acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1485e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1486e2446a98SMatthew Knepley    no intervening spaces.
1487e2446a98SMatthew Knepley 
14883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1489e2446a98SMatthew Knepley 
1490e2446a98SMatthew Knepley    Input Parameters:
1491e2446a98SMatthew Knepley +  opt - the option one is seeking
1492e2446a98SMatthew Knepley .  text - short string describing option
1493e2446a98SMatthew Knepley .  man - manual page for option
1494e2446a98SMatthew Knepley -  nmax - maximum number of values
1495e2446a98SMatthew Knepley 
1496e2446a98SMatthew Knepley    Output Parameter:
1497e2446a98SMatthew Knepley +  value - location to copy values
1498e2446a98SMatthew Knepley .  nmax - actual number of values found
1499e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1500e2446a98SMatthew Knepley 
1501e2446a98SMatthew Knepley    Level: beginner
1502e2446a98SMatthew Knepley 
1503e2446a98SMatthew Knepley    Notes:
1504e2446a98SMatthew Knepley    The user should pass in an array of doubles
1505e2446a98SMatthew Knepley 
1506e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1507e2446a98SMatthew Knepley 
1508e2446a98SMatthew Knepley    Concepts: options database^array of strings
1509e2446a98SMatthew Knepley 
1510e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1511acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1512e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1513e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1514acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1515a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1516e2446a98SMatthew Knepley @*/
15178c34d3f5SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1518e2446a98SMatthew Knepley {
1519e2446a98SMatthew Knepley   PetscErrorCode ierr;
1520e2446a98SMatthew Knepley   PetscInt       i;
15218c34d3f5SBarry Smith   PetscOption    amsopt;
1522e2446a98SMatthew Knepley 
1523e2446a98SMatthew Knepley   PetscFunctionBegin;
1524e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1525ace3abfcSBarry Smith     PetscBool *vals;
15261ae3d29cSBarry Smith 
152783355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
15281a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1529ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
15301ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
15311ae3d29cSBarry Smith     amsopt->arraylength = *n;
15321ae3d29cSBarry Smith   }
1533e55864a3SBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1534e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1535e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1536e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1537e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1538e2446a98SMatthew Knepley     }
1539e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1540e2446a98SMatthew Knepley   }
1541e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1542e2446a98SMatthew Knepley }
1543e2446a98SMatthew Knepley 
15448cc676e6SMatthew G Knepley #undef __FUNCT__
1545e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
15468cc676e6SMatthew G Knepley /*@C
1547d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
15488cc676e6SMatthew G Knepley 
15498cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15508cc676e6SMatthew G Knepley 
15518cc676e6SMatthew G Knepley    Input Parameters:
15528cc676e6SMatthew G Knepley +  opt - option name
15538cc676e6SMatthew G Knepley .  text - short string that describes the option
15548cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15558cc676e6SMatthew G Knepley 
15568cc676e6SMatthew G Knepley    Output Parameter:
15578cc676e6SMatthew G Knepley +  viewer - the viewer
15588cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15598cc676e6SMatthew G Knepley 
15608cc676e6SMatthew G Knepley    Level: beginner
15618cc676e6SMatthew G Knepley 
15628cc676e6SMatthew G Knepley    Concepts: options database^has int
15638cc676e6SMatthew G Knepley 
15648cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15658cc676e6SMatthew G Knepley 
15665a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
15678cc676e6SMatthew G Knepley 
15688cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15698cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15708cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15718cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15728cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15738cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1574a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15758cc676e6SMatthew G Knepley @*/
15768c34d3f5SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptions *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15778cc676e6SMatthew G Knepley {
15788cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15798c34d3f5SBarry Smith   PetscOption    amsopt;
15808cc676e6SMatthew G Knepley 
15818cc676e6SMatthew G Knepley   PetscFunctionBegin;
15821a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
15831a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
158464facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15855b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15868cc676e6SMatthew G Knepley   }
1587e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1588e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1589e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15908cc676e6SMatthew G Knepley   }
15918cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15928cc676e6SMatthew G Knepley }
15938cc676e6SMatthew G Knepley 
159453acd3b1SBarry Smith 
159553acd3b1SBarry Smith #undef __FUNCT__
159653acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
159753acd3b1SBarry Smith /*@C
1598b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
159953acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
160053acd3b1SBarry Smith 
16013f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
160253acd3b1SBarry Smith 
160353acd3b1SBarry Smith    Input Parameter:
160453acd3b1SBarry Smith .   head - the heading text
160553acd3b1SBarry Smith 
160653acd3b1SBarry Smith 
160753acd3b1SBarry Smith    Level: intermediate
160853acd3b1SBarry Smith 
160953acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
161053acd3b1SBarry Smith 
1611b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
161253acd3b1SBarry Smith 
161353acd3b1SBarry Smith    Concepts: options database^subheading
161453acd3b1SBarry Smith 
161553acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1616acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
161753acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
161853acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1619acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1620a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
162153acd3b1SBarry Smith @*/
16228c34d3f5SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptions *PetscOptionsObject,const char head[])
162353acd3b1SBarry Smith {
162453acd3b1SBarry Smith   PetscErrorCode ierr;
162553acd3b1SBarry Smith 
162653acd3b1SBarry Smith   PetscFunctionBegin;
1627e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1628e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
162953acd3b1SBarry Smith   }
163053acd3b1SBarry Smith   PetscFunctionReturn(0);
163153acd3b1SBarry Smith }
163253acd3b1SBarry Smith 
163353acd3b1SBarry Smith 
163453acd3b1SBarry Smith 
163553acd3b1SBarry Smith 
163653acd3b1SBarry Smith 
163753acd3b1SBarry Smith 
1638