xref: /petsc/src/sys/objects/aoptions.c (revision 1a1499c8e13c12f02cf4c59cfd6b0cfcce01ae9b)
17d0a6c19SBarry Smith 
2*1a1499c8SBarry 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 */
27e55864a3SBarry Smith PetscErrorCode PetscOptionsBegin_Private(PetscOptionsObjectType *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 */
53e55864a3SBarry Smith PetscErrorCode PetscObjectOptionsBegin_Private(PetscOptionsObjectType *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"
79e55864a3SBarry Smith static int PetscOptionsCreate_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscOptionType t,PetscOptions *amsopt)
8053acd3b1SBarry Smith {
8153acd3b1SBarry Smith   int          ierr;
8253acd3b1SBarry Smith   PetscOptions 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 
147aca04275SShao-Ching Huang #undef __FUNCT__
148aca04275SShao-Ching Huang #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 */
189e55864a3SBarry Smith PetscErrorCode PetscOptionsGetFromTextInput(PetscOptionsObjectType *PetscOptionsObject)
1906356e834SBarry Smith {
1916356e834SBarry Smith   PetscErrorCode ierr;
192e55864a3SBarry Smith   PetscOptions   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 
199e55864a3SBarry Smith   ierr = (*PetscPrintf)(PETSC_COMM_WORLD,"%s -------------------------------------------------\n",PetscOptionsObject->title);CHKERRQ(ierr);
2006356e834SBarry Smith   while (next) {
2016356e834SBarry Smith     switch (next->type) {
2026356e834SBarry Smith     case OPTION_HEAD:
2036356e834SBarry Smith       break;
204e26ddf31SBarry Smith     case OPTION_INT_ARRAY:
205e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
206e26ddf31SBarry Smith       vald = (PetscInt*) next->data;
207e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
208e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%d",vald[i]);CHKERRQ(ierr);
209e26ddf31SBarry Smith         if (i < next->arraylength-1) {
210e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
211e26ddf31SBarry Smith         }
212e26ddf31SBarry Smith       }
213e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
214e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
215e26ddf31SBarry Smith       if (str[0]) {
216e26ddf31SBarry Smith         PetscToken token;
217e26ddf31SBarry Smith         PetscInt   n=0,nmax = next->arraylength,*dvalue = (PetscInt*)next->data,start,end;
218e26ddf31SBarry Smith         size_t     len;
219e26ddf31SBarry Smith         char       *value;
220ace3abfcSBarry Smith         PetscBool  foundrange;
221e26ddf31SBarry Smith 
222e26ddf31SBarry Smith         next->set = PETSC_TRUE;
223e26ddf31SBarry Smith         value     = str;
224e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
225e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
226e26ddf31SBarry Smith         while (n < nmax) {
227e26ddf31SBarry Smith           if (!value) break;
228e26ddf31SBarry Smith 
229e26ddf31SBarry Smith           /* look for form  d-D where d and D are integers */
230e26ddf31SBarry Smith           foundrange = PETSC_FALSE;
231e26ddf31SBarry Smith           ierr       = PetscStrlen(value,&len);CHKERRQ(ierr);
232e26ddf31SBarry Smith           if (value[0] == '-') i=2;
233e26ddf31SBarry Smith           else i=1;
234330cf3c9SBarry Smith           for (;i<len; i++) {
235e26ddf31SBarry Smith             if (value[i] == '-') {
236e32f2f54SBarry Smith               if (i == len-1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Error in %D-th array entry %s\n",n,value);
237e26ddf31SBarry Smith               value[i] = 0;
238cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value,&start);CHKERRQ(ierr);
239cfbddea1SSatish Balay               ierr     = PetscOptionsStringToInt(value+i+1,&end);CHKERRQ(ierr);
240e32f2f54SBarry 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);
241e32f2f54SBarry 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);
242e26ddf31SBarry Smith               for (; start<end; start++) {
243e26ddf31SBarry Smith                 *dvalue = start; dvalue++;n++;
244e26ddf31SBarry Smith               }
245e26ddf31SBarry Smith               foundrange = PETSC_TRUE;
246e26ddf31SBarry Smith               break;
247e26ddf31SBarry Smith             }
248e26ddf31SBarry Smith           }
249e26ddf31SBarry Smith           if (!foundrange) {
250cfbddea1SSatish Balay             ierr = PetscOptionsStringToInt(value,dvalue);CHKERRQ(ierr);
251e26ddf31SBarry Smith             dvalue++;
252e26ddf31SBarry Smith             n++;
253e26ddf31SBarry Smith           }
254e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
255e26ddf31SBarry Smith         }
2568c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
257e26ddf31SBarry Smith       }
258e26ddf31SBarry Smith       break;
259e26ddf31SBarry Smith     case OPTION_REAL_ARRAY:
260e55864a3SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"-%s%s <",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",next->option+1);CHKERRQ(ierr);
261e26ddf31SBarry Smith       valr = (PetscReal*) next->data;
262e26ddf31SBarry Smith       for (i=0; i<next->arraylength; i++) {
263e26ddf31SBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"%g",valr[i]);CHKERRQ(ierr);
264e26ddf31SBarry Smith         if (i < next->arraylength-1) {
265e26ddf31SBarry Smith           ierr = PetscPrintf(PETSC_COMM_WORLD,",");CHKERRQ(ierr);
266e26ddf31SBarry Smith         }
267e26ddf31SBarry Smith       }
268e26ddf31SBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,">: %s (%s) ",next->text,next->man);CHKERRQ(ierr);
269e26ddf31SBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
270e26ddf31SBarry Smith       if (str[0]) {
271e26ddf31SBarry Smith         PetscToken token;
272e26ddf31SBarry Smith         PetscInt   n = 0,nmax = next->arraylength;
273e26ddf31SBarry Smith         PetscReal  *dvalue = (PetscReal*)next->data;
274e26ddf31SBarry Smith         char       *value;
275e26ddf31SBarry Smith 
276e26ddf31SBarry Smith         next->set = PETSC_TRUE;
277e26ddf31SBarry Smith         value     = str;
278e26ddf31SBarry Smith         ierr      = PetscTokenCreate(value,',',&token);CHKERRQ(ierr);
279e26ddf31SBarry Smith         ierr      = PetscTokenFind(token,&value);CHKERRQ(ierr);
280e26ddf31SBarry Smith         while (n < nmax) {
281e26ddf31SBarry Smith           if (!value) break;
282cfbddea1SSatish Balay           ierr = PetscOptionsStringToReal(value,dvalue);CHKERRQ(ierr);
283e26ddf31SBarry Smith           dvalue++;
284e26ddf31SBarry Smith           n++;
285e26ddf31SBarry Smith           ierr = PetscTokenFind(token,&value);CHKERRQ(ierr);
286e26ddf31SBarry Smith         }
2878c74ee41SBarry Smith         ierr = PetscTokenDestroy(&token);CHKERRQ(ierr);
288e26ddf31SBarry Smith       }
289e26ddf31SBarry Smith       break;
2906356e834SBarry Smith     case OPTION_INT:
291e55864a3SBarry 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);
2923fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
2933fc1eb6aSBarry Smith       if (str[0]) {
294d25d7f95SJed Brown #if defined(PETSC_SIZEOF_LONG_LONG)
295d25d7f95SJed Brown         long long lid;
296d25d7f95SJed Brown         sscanf(str,"%lld",&lid);
297*1a1499c8SBarry 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);
298c272547aSJed Brown #else
299d25d7f95SJed Brown         long  lid;
300d25d7f95SJed Brown         sscanf(str,"%ld",&lid);
301*1a1499c8SBarry 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);
302c272547aSJed Brown #endif
303a297a907SKarl Rupp 
304d25d7f95SJed Brown         next->set = PETSC_TRUE;
305d25d7f95SJed Brown         *((PetscInt*)next->data) = (PetscInt)lid;
306aee2cecaSBarry Smith       }
307aee2cecaSBarry Smith       break;
308aee2cecaSBarry Smith     case OPTION_REAL:
309e55864a3SBarry 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);
3103fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3113fc1eb6aSBarry Smith       if (str[0]) {
312ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
313a4404d99SBarry Smith         sscanf(str,"%e",&ir);
314ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
315aee2cecaSBarry Smith         sscanf(str,"%le",&ir);
316ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
317d9822059SBarry Smith         ir = strtoflt128(str,0);
318d9822059SBarry Smith #else
319513dbe71SLisandro Dalcin         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unknown scalar type");
320a4404d99SBarry Smith #endif
321aee2cecaSBarry Smith         next->set                 = PETSC_TRUE;
322aee2cecaSBarry Smith         *((PetscReal*)next->data) = ir;
323aee2cecaSBarry Smith       }
324aee2cecaSBarry Smith       break;
3257781c08eSBarry Smith     case OPTION_BOOL:
32683355fc5SBarry 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);
3277781c08eSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3287781c08eSBarry Smith       if (str[0]) {
3297781c08eSBarry Smith         ierr = PetscOptionsStringToBool(str,&bid);CHKERRQ(ierr);
3307781c08eSBarry Smith         next->set = PETSC_TRUE;
3317781c08eSBarry Smith         *((PetscBool*)next->data) = bid;
3327781c08eSBarry Smith       }
3337781c08eSBarry Smith       break;
334aee2cecaSBarry Smith     case OPTION_STRING:
335e55864a3SBarry 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);
3363fc1eb6aSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3373fc1eb6aSBarry Smith       if (str[0]) {
338aee2cecaSBarry Smith         next->set = PETSC_TRUE;
33964facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3405b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3416356e834SBarry Smith       }
3426356e834SBarry Smith       break;
343a264d7a6SBarry Smith     case OPTION_FLIST:
344e55864a3SBarry Smith       ierr = PetscFunctionListPrintTypes(PETSC_COMM_WORLD,stdout,PetscOptionsObject->prefix,next->option,next->text,next->man,next->flist,(char*)next->data);CHKERRQ(ierr);
3453cc1e11dSBarry Smith       ierr = PetscScanString(PETSC_COMM_WORLD,512,str);CHKERRQ(ierr);
3463cc1e11dSBarry Smith       if (str[0]) {
347e55864a3SBarry Smith         PetscOptionsObject->changedmethod = PETSC_TRUE;
3483cc1e11dSBarry Smith         next->set = PETSC_TRUE;
34964facd6cSBarry Smith         /* must use system malloc since SAWs may free this */
3505b02f95dSBarry Smith         ierr = PetscStrdup(str,(char**)&next->data);CHKERRQ(ierr);
3513cc1e11dSBarry Smith       }
3523cc1e11dSBarry Smith       break;
353b432afdaSMatthew Knepley     default:
354b432afdaSMatthew Knepley       break;
3556356e834SBarry Smith     }
3566356e834SBarry Smith     next = next->next;
3576356e834SBarry Smith   }
3586356e834SBarry Smith   PetscFunctionReturn(0);
3596356e834SBarry Smith }
3606356e834SBarry Smith 
361e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
362e04113cfSBarry Smith #include <petscviewersaws.h>
363d5649816SBarry Smith 
364d5649816SBarry Smith static int count = 0;
365d5649816SBarry Smith 
366b3506946SBarry Smith #undef __FUNCT__
367e04113cfSBarry Smith #define __FUNCT__ "PetscOptionsSAWsDestroy"
368e04113cfSBarry Smith PetscErrorCode PetscOptionsSAWsDestroy(void)
369d5649816SBarry Smith {
3702657e9d9SBarry Smith   PetscFunctionBegin;
371d5649816SBarry Smith   PetscFunctionReturn(0);
372d5649816SBarry Smith }
373d5649816SBarry Smith 
374d5649816SBarry Smith #undef __FUNCT__
3757781c08eSBarry Smith #define __FUNCT__ "PetscOptionsSAWsInput"
376b3506946SBarry Smith /*
3777781c08eSBarry Smith     PetscOptionsSAWsInput - Presents all the PETSc Options processed by the program so the user may change them at runtime using the SAWs
378b3506946SBarry Smith 
379b3506946SBarry Smith     Bugs:
380b3506946SBarry Smith +    All processes must traverse through the exact same set of option queries do to the call to PetscScanString()
381b3506946SBarry Smith .    Internal strings have arbitrary length and string copies are not checked that they fit into string space
382b3506946SBarry Smith -    Only works for PetscInt == int, PetscReal == double etc
383b3506946SBarry Smith 
384b3506946SBarry Smith 
385b3506946SBarry Smith */
386e55864a3SBarry Smith PetscErrorCode PetscOptionsAMSInput(PetscOptionsObjectType *PetscOptionsObject)
387b3506946SBarry Smith {
388b3506946SBarry Smith   PetscErrorCode ierr;
389e55864a3SBarry Smith   PetscOptions   next     = PetscOptionsObject->next;
390d5649816SBarry Smith   static int     mancount = 0;
391b3506946SBarry Smith   char           options[16];
3927aab2a10SBarry Smith   PetscBool      changedmethod = PETSC_FALSE;
39388a9752cSBarry Smith   char           manname[16],textname[16];
3942657e9d9SBarry Smith   char           dir[1024];
395b3506946SBarry Smith 
396b3506946SBarry Smith   /* the next line is a bug, this will only work if all processors are here, the comm passed in is ignored!!! */
397b3506946SBarry Smith   sprintf(options,"Options_%d",count++);
398a297a907SKarl Rupp 
399e55864a3SBarry Smith   PetscOptionsObject->pprefix = PetscOptionsObject->prefix; /* AMS will change this, so cannot pass prefix directly */
4001bc75a8dSBarry Smith 
401d50831c4SBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","_title");CHKERRQ(ierr);
40283355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->title,1,SAWs_READ,SAWs_STRING));
4037781c08eSBarry Smith   ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s","prefix");CHKERRQ(ierr);
40483355fc5SBarry Smith   PetscStackCallSAWs(SAWs_Register,(dir,&PetscOptionsObject->pprefix,1,SAWs_READ,SAWs_STRING));
4052657e9d9SBarry Smith   PetscStackCallSAWs(SAWs_Register,("/PETSc/Options/ChangedMethod",&changedmethod,1,SAWs_WRITE,SAWs_BOOLEAN));
406b3506946SBarry Smith 
407b3506946SBarry Smith   while (next) {
408d50831c4SBarry Smith     sprintf(manname,"_man_%d",mancount);
4092657e9d9SBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",manname);CHKERRQ(ierr);
4107781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->man,1,SAWs_READ,SAWs_STRING));
411d50831c4SBarry Smith     sprintf(textname,"_text_%d",mancount++);
4127781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",textname);CHKERRQ(ierr);
4137781c08eSBarry Smith     PetscStackCallSAWs(SAWs_Register,(dir,&next->text,1,SAWs_READ,SAWs_STRING));
4149f32e415SBarry Smith 
415b3506946SBarry Smith     switch (next->type) {
416b3506946SBarry Smith     case OPTION_HEAD:
417b3506946SBarry Smith       break;
418b3506946SBarry Smith     case OPTION_INT_ARRAY:
4197781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4202657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_INT));
421b3506946SBarry Smith       break;
422b3506946SBarry Smith     case OPTION_REAL_ARRAY:
4237781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4242657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_DOUBLE));
425b3506946SBarry Smith       break;
426b3506946SBarry Smith     case OPTION_INT:
4277781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4282657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_INT));
429b3506946SBarry Smith       break;
430b3506946SBarry Smith     case OPTION_REAL:
4317781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4322657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_DOUBLE));
433b3506946SBarry Smith       break;
4347781c08eSBarry Smith     case OPTION_BOOL:
4357781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4362657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,1,SAWs_WRITE,SAWs_BOOLEAN));
4371ae3d29cSBarry Smith       break;
4387781c08eSBarry Smith     case OPTION_BOOL_ARRAY:
4397781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4402657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_BOOLEAN));
44171f08665SBarry Smith       break;
442b3506946SBarry Smith     case OPTION_STRING:
4437781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4447781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
4451ae3d29cSBarry Smith       break;
4461ae3d29cSBarry Smith     case OPTION_STRING_ARRAY:
4477781c08eSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4482657e9d9SBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,next->data,next->arraylength,SAWs_WRITE,SAWs_STRING));
449b3506946SBarry Smith       break;
450a264d7a6SBarry Smith     case OPTION_FLIST:
451a264d7a6SBarry Smith       {
452a264d7a6SBarry Smith       PetscInt ntext;
4537781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4547781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
455a264d7a6SBarry Smith       ierr = PetscFunctionListGet(next->flist,(const char***)&next->edata,&ntext);CHKERRQ(ierr);
456a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
457a264d7a6SBarry Smith       }
458a264d7a6SBarry Smith       break;
4591ae3d29cSBarry Smith     case OPTION_ELIST:
460a264d7a6SBarry Smith       {
461a264d7a6SBarry Smith       PetscInt ntext = next->nlist;
4627781c08eSBarry Smith       ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
4637781c08eSBarry Smith       PetscStackCallSAWs(SAWs_Register,(dir,&next->data,1,SAWs_WRITE,SAWs_STRING));
464315e5ce4SBarry Smith       ierr = PetscMalloc((ntext+1)*sizeof(char),&next->edata);CHKERRQ(ierr);
4651ae3d29cSBarry Smith       ierr = PetscMemcpy(next->edata,next->list,ntext*sizeof(char*));CHKERRQ(ierr);
466a264d7a6SBarry Smith       PetscStackCallSAWs(SAWs_Set_Legal_Variable_Values,(dir,ntext,next->edata));
467a264d7a6SBarry Smith       }
468a264d7a6SBarry Smith       break;
469b3506946SBarry Smith     default:
470b3506946SBarry Smith       break;
471b3506946SBarry Smith     }
472b3506946SBarry Smith     next = next->next;
473b3506946SBarry Smith   }
474b3506946SBarry Smith 
475b3506946SBarry Smith   /* wait until accessor has unlocked the memory */
4767aab2a10SBarry Smith   ierr = PetscSAWsBlock();CHKERRQ(ierr);
477b3506946SBarry Smith 
47888a9752cSBarry Smith   /* determine if any values have been set in GUI */
47983355fc5SBarry Smith   next = PetscOptionsObject->next;
48088a9752cSBarry Smith   while (next) {
48188a9752cSBarry Smith     ierr = PetscSNPrintf(dir,1024,"/PETSc/Options/%s",next->option);CHKERRQ(ierr);
48288a9752cSBarry Smith     PetscStackCallSAWs(SAWs_Selected,(dir,&next->set));
48388a9752cSBarry Smith     next = next->next;
48488a9752cSBarry Smith   }
48588a9752cSBarry Smith 
486b3506946SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
487b3506946SBarry Smith   if (changedmethod) PetscOptionsPublishCount = -2;
488b3506946SBarry Smith 
4899a492a5cSBarry Smith   PetscStackCallSAWs(SAWs_Delete,("/PETSc/Options"));
490b3506946SBarry Smith   PetscFunctionReturn(0);
491b3506946SBarry Smith }
492b3506946SBarry Smith #endif
493b3506946SBarry Smith 
4946356e834SBarry Smith #undef __FUNCT__
49553acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsEnd_Private"
496e55864a3SBarry Smith PetscErrorCode PetscOptionsEnd_Private(PetscOptionsObjectType *PetscOptionsObject)
49753acd3b1SBarry Smith {
49853acd3b1SBarry Smith   PetscErrorCode ierr;
4996356e834SBarry Smith   PetscOptions   last;
5006356e834SBarry Smith   char           option[256],value[1024],tmp[32];
501330cf3c9SBarry Smith   size_t         j;
50253acd3b1SBarry Smith 
50353acd3b1SBarry Smith   PetscFunctionBegin;
50483355fc5SBarry Smith   if (PetscOptionsObject->next) {
50583355fc5SBarry Smith     if (!PetscOptionsObject->count) {
506a264d7a6SBarry Smith #if defined(PETSC_HAVE_SAWS)
507475446a1SBarry Smith       ierr = PetscOptionsSAWsInput();CHKERRQ(ierr);
508b3506946SBarry Smith #else
509e55864a3SBarry Smith       ierr = PetscOptionsGetFromTextInput(PetscOptionsObject);CHKERRQ(ierr);
510b3506946SBarry Smith #endif
511aee2cecaSBarry Smith     }
512aee2cecaSBarry Smith   }
5136356e834SBarry Smith 
514e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->title);CHKERRQ(ierr);
515e55864a3SBarry Smith   ierr = PetscFree(PetscOptionsObject->prefix);CHKERRQ(ierr);
5166356e834SBarry Smith 
5176356e834SBarry Smith   /* reset counter to -2; this updates the screen with the new options for the selected method */
518e55864a3SBarry Smith   if (PetscOptionsObject->changedmethod) PetscOptionsObject->count = -2;
51961b37b28SSatish Balay   /* reset alreadyprinted flag */
520e55864a3SBarry Smith   PetscOptionsObject->alreadyprinted = PETSC_FALSE;
521e55864a3SBarry Smith   if (PetscOptionsObject->object) PetscOptionsObject->object->optionsprinted = PETSC_TRUE;
522e55864a3SBarry Smith   PetscOptionsObject->object = NULL;
5236356e834SBarry Smith 
524e55864a3SBarry Smith   while (PetscOptionsObject->next) {
525e55864a3SBarry Smith     if (PetscOptionsObject->next->set) {
526e55864a3SBarry Smith       if (PetscOptionsObject->prefix) {
5276356e834SBarry Smith         ierr = PetscStrcpy(option,"-");CHKERRQ(ierr);
528e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->prefix);CHKERRQ(ierr);
529e55864a3SBarry Smith         ierr = PetscStrcat(option,PetscOptionsObject->next->option+1);CHKERRQ(ierr);
5306356e834SBarry Smith       } else {
531e55864a3SBarry Smith         ierr = PetscStrcpy(option,PetscOptionsObject->next->option);CHKERRQ(ierr);
5326356e834SBarry Smith       }
5336356e834SBarry Smith 
534e55864a3SBarry Smith       switch (PetscOptionsObject->next->type) {
5356356e834SBarry Smith       case OPTION_HEAD:
5366356e834SBarry Smith         break;
537e26ddf31SBarry Smith       case OPTION_INT_ARRAY:
538e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[0]);
539e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
540e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscInt*)PetscOptionsObject->next->data)[j]);
541e26ddf31SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
542e26ddf31SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
543e26ddf31SBarry Smith         }
544e26ddf31SBarry Smith         break;
5456356e834SBarry Smith       case OPTION_INT:
546e55864a3SBarry Smith         sprintf(value,"%d",(int) *(PetscInt*)PetscOptionsObject->next->data);
5476356e834SBarry Smith         break;
5486356e834SBarry Smith       case OPTION_REAL:
549e55864a3SBarry Smith         sprintf(value,"%g",(double) *(PetscReal*)PetscOptionsObject->next->data);
5506356e834SBarry Smith         break;
5516356e834SBarry Smith       case OPTION_REAL_ARRAY:
552e55864a3SBarry Smith         sprintf(value,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[0]);
553e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
554e55864a3SBarry Smith           sprintf(tmp,"%g",(double)((PetscReal*)PetscOptionsObject->next->data)[j]);
5556356e834SBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5566356e834SBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5576356e834SBarry Smith         }
5586356e834SBarry Smith         break;
5597781c08eSBarry Smith       case OPTION_BOOL:
560e55864a3SBarry Smith         sprintf(value,"%d",*(int*)PetscOptionsObject->next->data);
5616356e834SBarry Smith         break;
5627781c08eSBarry Smith       case OPTION_BOOL_ARRAY:
563e55864a3SBarry Smith         sprintf(value,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[0]);
564e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
565e55864a3SBarry Smith           sprintf(tmp,"%d",(int)((PetscBool*)PetscOptionsObject->next->data)[j]);
5661ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5671ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5681ae3d29cSBarry Smith         }
5691ae3d29cSBarry Smith         break;
570a264d7a6SBarry Smith       case OPTION_FLIST:
5711ae3d29cSBarry Smith       case OPTION_ELIST:
572e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
5736356e834SBarry Smith         break;
5741ae3d29cSBarry Smith       case OPTION_STRING:
575e55864a3SBarry Smith         ierr = PetscStrcpy(value,(char*)PetscOptionsObject->next->data);CHKERRQ(ierr);
5761ae3d29cSBarry Smith       case OPTION_STRING_ARRAY:
577e55864a3SBarry Smith         sprintf(value,"%s",((char**)PetscOptionsObject->next->data)[0]);
578e55864a3SBarry Smith         for (j=1; j<PetscOptionsObject->next->arraylength; j++) {
579e55864a3SBarry Smith           sprintf(tmp,"%s",((char**)PetscOptionsObject->next->data)[j]);
5801ae3d29cSBarry Smith           ierr = PetscStrcat(value,",");CHKERRQ(ierr);
5811ae3d29cSBarry Smith           ierr = PetscStrcat(value,tmp);CHKERRQ(ierr);
5821ae3d29cSBarry Smith         }
5836356e834SBarry Smith         break;
5846356e834SBarry Smith       }
5856356e834SBarry Smith       ierr = PetscOptionsSetValue(option,value);CHKERRQ(ierr);
5866356e834SBarry Smith     }
587e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->text);CHKERRQ(ierr);
588e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->option);CHKERRQ(ierr);
589e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->man);CHKERRQ(ierr);
590e55864a3SBarry Smith     ierr   = PetscFree(PetscOptionsObject->next->edata);CHKERRQ(ierr);
591c979a496SBarry Smith 
59283355fc5SBarry Smith     if ((PetscOptionsObject->next->type == OPTION_STRING) || (PetscOptionsObject->next->type == OPTION_FLIST) || (PetscOptionsObject->next->type == OPTION_ELIST)){
59383355fc5SBarry Smith       free(PetscOptionsObject->next->data);
594c979a496SBarry Smith     } else {
59583355fc5SBarry Smith       ierr   = PetscFree(PetscOptionsObject->next->data);CHKERRQ(ierr);
596c979a496SBarry Smith     }
5977781c08eSBarry Smith 
59883355fc5SBarry Smith     last                    = PetscOptionsObject->next;
59983355fc5SBarry Smith     PetscOptionsObject->next = PetscOptionsObject->next->next;
6006356e834SBarry Smith     ierr                    = PetscFree(last);CHKERRQ(ierr);
6016356e834SBarry Smith   }
602e55864a3SBarry Smith   PetscOptionsObject->next = 0;
60353acd3b1SBarry Smith   PetscFunctionReturn(0);
60453acd3b1SBarry Smith }
60553acd3b1SBarry Smith 
60653acd3b1SBarry Smith #undef __FUNCT__
607e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEnum_Private"
60853acd3b1SBarry Smith /*@C
60953acd3b1SBarry Smith    PetscOptionsEnum - Gets the enum value for a particular option in the database.
61053acd3b1SBarry Smith 
6113f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
61253acd3b1SBarry Smith 
61353acd3b1SBarry Smith    Input Parameters:
61453acd3b1SBarry Smith +  opt - option name
61553acd3b1SBarry Smith .  text - short string that describes the option
61653acd3b1SBarry Smith .  man - manual page with additional information on option
61753acd3b1SBarry Smith .  list - array containing the list of choices, followed by the enum name, followed by the enum prefix, followed by a null
6180fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6190fdccdaeSBarry Smith $                 PetscOptionsEnum(..., obj->value,&object->value,...) or
6200fdccdaeSBarry Smith $                 value = defaultvalue
6210fdccdaeSBarry Smith $                 PetscOptionsEnum(..., value,&value,&flg);
6220fdccdaeSBarry Smith $                 if (flg) {
62353acd3b1SBarry Smith 
62453acd3b1SBarry Smith    Output Parameter:
62553acd3b1SBarry Smith +  value - the  value to return
626b32e0204SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
62753acd3b1SBarry Smith 
62853acd3b1SBarry Smith    Level: beginner
62953acd3b1SBarry Smith 
63053acd3b1SBarry Smith    Concepts: options database
63153acd3b1SBarry Smith 
63253acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
63353acd3b1SBarry Smith 
63453acd3b1SBarry Smith           list is usually something like PCASMTypes or some other predefined list of enum names
63553acd3b1SBarry Smith 
63653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
637acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
638acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
63953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
64053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
641acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
642a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
64353acd3b1SBarry Smith @*/
644*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsEnum_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],const char *const *list,PetscEnum currentvalue,PetscEnum *value,PetscBool  *set)
64553acd3b1SBarry Smith {
64653acd3b1SBarry Smith   PetscErrorCode ierr;
64753acd3b1SBarry Smith   PetscInt       ntext = 0;
648aa5bb8c0SSatish Balay   PetscInt       tval;
649ace3abfcSBarry Smith   PetscBool      tflg;
65053acd3b1SBarry Smith 
65153acd3b1SBarry Smith   PetscFunctionBegin;
65253acd3b1SBarry Smith   while (list[ntext++]) {
653e32f2f54SBarry Smith     if (ntext > 50) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument appears to be wrong or have more than 50 entries");
65453acd3b1SBarry Smith   }
655e32f2f54SBarry Smith   if (ntext < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"List argument must have at least two entries: typename and type prefix");
65653acd3b1SBarry Smith   ntext -= 3;
657e55864a3SBarry Smith   ierr   = PetscOptionsEList_Private(PetscOptionsObject,opt,text,man,list,ntext,list[currentvalue],&tval,&tflg);CHKERRQ(ierr);
658aa5bb8c0SSatish Balay   /* with PETSC_USE_64BIT_INDICES sizeof(PetscInt) != sizeof(PetscEnum) */
659aa5bb8c0SSatish Balay   if (tflg) *value = (PetscEnum)tval;
660aa5bb8c0SSatish Balay   if (set)  *set   = tflg;
66153acd3b1SBarry Smith   PetscFunctionReturn(0);
66253acd3b1SBarry Smith }
66353acd3b1SBarry Smith 
66453acd3b1SBarry Smith /* -------------------------------------------------------------------------------------------------------------*/
66553acd3b1SBarry Smith #undef __FUNCT__
666e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsInt_Private"
66753acd3b1SBarry Smith /*@C
66853acd3b1SBarry Smith    PetscOptionsInt - Gets the integer value for a particular option in the database.
66953acd3b1SBarry Smith 
6703f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
67153acd3b1SBarry Smith 
67253acd3b1SBarry Smith    Input Parameters:
67353acd3b1SBarry Smith +  opt - option name
67453acd3b1SBarry Smith .  text - short string that describes the option
67553acd3b1SBarry Smith .  man - manual page with additional information on option
6760fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
6770fdccdaeSBarry Smith $                 PetscOptionsInt(..., obj->value,&object->value,...) or
6780fdccdaeSBarry Smith $                 value = defaultvalue
6790fdccdaeSBarry Smith $                 PetscOptionsInt(..., value,&value,&flg);
6800fdccdaeSBarry Smith $                 if (flg) {
68153acd3b1SBarry Smith 
68253acd3b1SBarry Smith    Output Parameter:
68353acd3b1SBarry Smith +  value - the integer value to return
68453acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
68553acd3b1SBarry Smith 
68653acd3b1SBarry Smith    Level: beginner
68753acd3b1SBarry Smith 
68853acd3b1SBarry Smith    Concepts: options database^has int
68953acd3b1SBarry Smith 
69053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
69153acd3b1SBarry Smith 
69253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
693acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
694acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
69553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
69653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
697acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
698a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
69953acd3b1SBarry Smith @*/
700*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsInt_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt currentvalue,PetscInt *value,PetscBool  *set)
70153acd3b1SBarry Smith {
70253acd3b1SBarry Smith   PetscErrorCode ierr;
7036356e834SBarry Smith   PetscOptions   amsopt;
70453acd3b1SBarry Smith 
70553acd3b1SBarry Smith   PetscFunctionBegin;
706e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
707e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT,&amsopt);CHKERRQ(ierr);
7086356e834SBarry Smith     ierr = PetscMalloc(sizeof(PetscInt),&amsopt->data);CHKERRQ(ierr);
709a297a907SKarl Rupp 
7100fdccdaeSBarry Smith     *(PetscInt*)amsopt->data = currentvalue;
711af6d86caSBarry Smith   }
712e55864a3SBarry Smith   ierr = PetscOptionsGetInt(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
713e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
714*1a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
71553acd3b1SBarry Smith   }
71653acd3b1SBarry Smith   PetscFunctionReturn(0);
71753acd3b1SBarry Smith }
71853acd3b1SBarry Smith 
71953acd3b1SBarry Smith #undef __FUNCT__
720e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsString_Private"
72153acd3b1SBarry Smith /*@C
72253acd3b1SBarry Smith    PetscOptionsString - Gets the string value for a particular option in the database.
72353acd3b1SBarry Smith 
7243f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
72553acd3b1SBarry Smith 
72653acd3b1SBarry Smith    Input Parameters:
72753acd3b1SBarry Smith +  opt - option name
72853acd3b1SBarry Smith .  text - short string that describes the option
72953acd3b1SBarry Smith .  man - manual page with additional information on option
7300fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. This is not used to set value
731bcbf2dc5SJed Brown -  len - length of the result string including null terminator
73253acd3b1SBarry Smith 
73353acd3b1SBarry Smith    Output Parameter:
73453acd3b1SBarry Smith +  value - the value to return
73553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
73653acd3b1SBarry Smith 
73753acd3b1SBarry Smith    Level: beginner
73853acd3b1SBarry Smith 
73953acd3b1SBarry Smith    Concepts: options database^has int
74053acd3b1SBarry Smith 
74153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
74253acd3b1SBarry Smith 
7437fccdfe4SBarry 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).
7447fccdfe4SBarry Smith 
74553acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
746acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
747acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsReal(), PetscOptionsBool(),
74853acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
74953acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
750acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
751a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
75253acd3b1SBarry Smith @*/
753*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsString_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],const char currentvalue[],char value[],size_t len,PetscBool  *set)
75453acd3b1SBarry Smith {
75553acd3b1SBarry Smith   PetscErrorCode ierr;
756aee2cecaSBarry Smith   PetscOptions   amsopt;
75753acd3b1SBarry Smith 
75853acd3b1SBarry Smith   PetscFunctionBegin;
759*1a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
760*1a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
76164facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
7620fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
763af6d86caSBarry Smith   }
764e55864a3SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
765e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
766*1a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,currentvalue,text,ManSection(man));CHKERRQ(ierr);
76753acd3b1SBarry Smith   }
76853acd3b1SBarry Smith   PetscFunctionReturn(0);
76953acd3b1SBarry Smith }
77053acd3b1SBarry Smith 
77153acd3b1SBarry Smith #undef __FUNCT__
772e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsReal_Private"
77353acd3b1SBarry Smith /*@C
77453acd3b1SBarry Smith    PetscOptionsReal - Gets the PetscReal value for a particular option in the database.
77553acd3b1SBarry Smith 
7763f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
77753acd3b1SBarry Smith 
77853acd3b1SBarry Smith    Input Parameters:
77953acd3b1SBarry Smith +  opt - option name
78053acd3b1SBarry Smith .  text - short string that describes the option
78153acd3b1SBarry Smith .  man - manual page with additional information on option
7820fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
7830fdccdaeSBarry Smith $                 PetscOptionsReal(..., obj->value,&object->value,...) or
7840fdccdaeSBarry Smith $                 value = defaultvalue
7850fdccdaeSBarry Smith $                 PetscOptionsReal(..., value,&value,&flg);
7860fdccdaeSBarry Smith $                 if (flg) {
78753acd3b1SBarry Smith 
78853acd3b1SBarry Smith    Output Parameter:
78953acd3b1SBarry Smith +  value - the value to return
79053acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
79153acd3b1SBarry Smith 
79253acd3b1SBarry Smith    Level: beginner
79353acd3b1SBarry Smith 
79453acd3b1SBarry Smith    Concepts: options database^has int
79553acd3b1SBarry Smith 
79653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
79753acd3b1SBarry Smith 
79853acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
799acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
800acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
80153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
80253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
803acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
804a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
80553acd3b1SBarry Smith @*/
806*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsReal_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal currentvalue,PetscReal *value,PetscBool  *set)
80753acd3b1SBarry Smith {
80853acd3b1SBarry Smith   PetscErrorCode ierr;
809538aa990SBarry Smith   PetscOptions   amsopt;
81053acd3b1SBarry Smith 
81153acd3b1SBarry Smith   PetscFunctionBegin;
812e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
813e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL,&amsopt);CHKERRQ(ierr);
814538aa990SBarry Smith     ierr = PetscMalloc(sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
815a297a907SKarl Rupp 
8160fdccdaeSBarry Smith     *(PetscReal*)amsopt->data = currentvalue;
817538aa990SBarry Smith   }
818*1a1499c8SBarry Smith   ierr = PetscOptionsGetReal(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
819*1a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
820*1a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%g>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,(double)currentvalue,text,ManSection(man));CHKERRQ(ierr);
82153acd3b1SBarry Smith   }
82253acd3b1SBarry Smith   PetscFunctionReturn(0);
82353acd3b1SBarry Smith }
82453acd3b1SBarry Smith 
82553acd3b1SBarry Smith #undef __FUNCT__
826e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsScalar_Private"
82753acd3b1SBarry Smith /*@C
82853acd3b1SBarry Smith    PetscOptionsScalar - Gets the scalar value for a particular option in the database.
82953acd3b1SBarry Smith 
8303f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
83153acd3b1SBarry Smith 
83253acd3b1SBarry Smith    Input Parameters:
83353acd3b1SBarry Smith +  opt - option name
83453acd3b1SBarry Smith .  text - short string that describes the option
83553acd3b1SBarry Smith .  man - manual page with additional information on option
8360fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with either
8370fdccdaeSBarry Smith $                 PetscOptionsScalar(..., obj->value,&object->value,...) or
8380fdccdaeSBarry Smith $                 value = defaultvalue
8390fdccdaeSBarry Smith $                 PetscOptionsScalar(..., value,&value,&flg);
8400fdccdaeSBarry Smith $                 if (flg) {
8410fdccdaeSBarry Smith 
84253acd3b1SBarry Smith 
84353acd3b1SBarry Smith    Output Parameter:
84453acd3b1SBarry Smith +  value - the value to return
84553acd3b1SBarry Smith -  flg - PETSC_TRUE if found, else PETSC_FALSE
84653acd3b1SBarry Smith 
84753acd3b1SBarry Smith    Level: beginner
84853acd3b1SBarry Smith 
84953acd3b1SBarry Smith    Concepts: options database^has int
85053acd3b1SBarry Smith 
85153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
85253acd3b1SBarry Smith 
85353acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
854acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
855acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
85653acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
85753acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
858acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
859a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
86053acd3b1SBarry Smith @*/
861*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsScalar_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscScalar currentvalue,PetscScalar *value,PetscBool  *set)
86253acd3b1SBarry Smith {
86353acd3b1SBarry Smith   PetscErrorCode ierr;
86453acd3b1SBarry Smith 
86553acd3b1SBarry Smith   PetscFunctionBegin;
86653acd3b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8670fdccdaeSBarry Smith   ierr = PetscOptionsReal(opt,text,man,currentvalue,value,set);CHKERRQ(ierr);
86853acd3b1SBarry Smith #else
869e55864a3SBarry Smith   ierr = PetscOptionsGetScalar(PetscOptionsObject->prefix,opt,value,set);CHKERRQ(ierr);
87053acd3b1SBarry Smith #endif
87153acd3b1SBarry Smith   PetscFunctionReturn(0);
87253acd3b1SBarry Smith }
87353acd3b1SBarry Smith 
87453acd3b1SBarry Smith #undef __FUNCT__
875e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsName_Private"
87653acd3b1SBarry Smith /*@C
87790d69ab7SBarry 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
87890d69ab7SBarry Smith                       its value is set to false.
87953acd3b1SBarry Smith 
8803f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
88153acd3b1SBarry Smith 
88253acd3b1SBarry Smith    Input Parameters:
88353acd3b1SBarry Smith +  opt - option name
88453acd3b1SBarry Smith .  text - short string that describes the option
88553acd3b1SBarry Smith -  man - manual page with additional information on option
88653acd3b1SBarry Smith 
88753acd3b1SBarry Smith    Output Parameter:
88853acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
88953acd3b1SBarry Smith 
89053acd3b1SBarry Smith    Level: beginner
89153acd3b1SBarry Smith 
89253acd3b1SBarry Smith    Concepts: options database^has int
89353acd3b1SBarry Smith 
89453acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
89553acd3b1SBarry Smith 
89653acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
897acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
898acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
89953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
90053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
901acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
902a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
90353acd3b1SBarry Smith @*/
904e55864a3SBarry Smith PetscErrorCode  PetscOptionsName_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
90553acd3b1SBarry Smith {
90653acd3b1SBarry Smith   PetscErrorCode ierr;
9071ae3d29cSBarry Smith   PetscOptions   amsopt;
90853acd3b1SBarry Smith 
90953acd3b1SBarry Smith   PetscFunctionBegin;
910e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
911e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
912ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
913a297a907SKarl Rupp 
914ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
9151ae3d29cSBarry Smith   }
916e55864a3SBarry Smith   ierr = PetscOptionsHasName(PetscOptionsObject->prefix,opt,flg);CHKERRQ(ierr);
917e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
918e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
91953acd3b1SBarry Smith   }
92053acd3b1SBarry Smith   PetscFunctionReturn(0);
92153acd3b1SBarry Smith }
92253acd3b1SBarry Smith 
92353acd3b1SBarry Smith #undef __FUNCT__
924e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsFList_Private"
92553acd3b1SBarry Smith /*@C
926a264d7a6SBarry Smith      PetscOptionsFList - Puts a list of option values that a single one may be selected from
92753acd3b1SBarry Smith 
9283f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
92953acd3b1SBarry Smith 
93053acd3b1SBarry Smith    Input Parameters:
93153acd3b1SBarry Smith +  opt - option name
93253acd3b1SBarry Smith .  text - short string that describes the option
93353acd3b1SBarry Smith .  man - manual page with additional information on option
93453acd3b1SBarry Smith .  list - the possible choices
9350fdccdaeSBarry Smith .  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
9360fdccdaeSBarry Smith $                 PetscOptionsFlist(..., obj->value,value,len,&flg);
9370fdccdaeSBarry Smith $                 if (flg) {
9383cc1e11dSBarry Smith -  len - the length of the character array value
93953acd3b1SBarry Smith 
94053acd3b1SBarry Smith    Output Parameter:
94153acd3b1SBarry Smith +  value - the value to return
94253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
94353acd3b1SBarry Smith 
94453acd3b1SBarry Smith    Level: intermediate
94553acd3b1SBarry Smith 
94653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
94753acd3b1SBarry Smith 
94853acd3b1SBarry Smith    See PetscOptionsEList() for when the choices are given in a string array
94953acd3b1SBarry Smith 
95053acd3b1SBarry Smith    To get a listing of all currently specified options,
95188c29154SBarry Smith     see PetscOptionsView() or PetscOptionsGetAll()
95253acd3b1SBarry Smith 
953eabe10d7SBarry Smith    Developer Note: This cannot check for invalid selection because of things like MATAIJ that are not included in the list
954eabe10d7SBarry Smith 
95553acd3b1SBarry Smith    Concepts: options database^list
95653acd3b1SBarry Smith 
95753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
958acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
95953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
96053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
961acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
962a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsEnum()
96353acd3b1SBarry Smith @*/
964*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsFList_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char ltext[],const char man[],PetscFunctionList list,const char currentvalue[],char value[],size_t len,PetscBool  *set)
96553acd3b1SBarry Smith {
96653acd3b1SBarry Smith   PetscErrorCode ierr;
9673cc1e11dSBarry Smith   PetscOptions   amsopt;
96853acd3b1SBarry Smith 
96953acd3b1SBarry Smith   PetscFunctionBegin;
970*1a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
971*1a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_FLIST,&amsopt);CHKERRQ(ierr);
97264facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
9730fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
9743cc1e11dSBarry Smith     amsopt->flist = list;
9753cc1e11dSBarry Smith   }
976*1a1499c8SBarry Smith   ierr = PetscOptionsGetString(PetscOptionsObject->prefix,opt,value,len,set);CHKERRQ(ierr);
977*1a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
978*1a1499c8SBarry Smith     ierr = PetscFunctionListPrintTypes(PetscOptionsObject->comm,stdout,PetscOptionsObject->prefix,opt,ltext,man,list,currentvalue);CHKERRQ(ierr);CHKERRQ(ierr);
97953acd3b1SBarry Smith   }
98053acd3b1SBarry Smith   PetscFunctionReturn(0);
98153acd3b1SBarry Smith }
98253acd3b1SBarry Smith 
98353acd3b1SBarry Smith #undef __FUNCT__
984e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsEList_Private"
98553acd3b1SBarry Smith /*@C
98653acd3b1SBarry Smith      PetscOptionsEList - Puts a list of option values that a single one may be selected from
98753acd3b1SBarry Smith 
9883f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
98953acd3b1SBarry Smith 
99053acd3b1SBarry Smith    Input Parameters:
99153acd3b1SBarry Smith +  opt - option name
99253acd3b1SBarry Smith .  ltext - short string that describes the option
99353acd3b1SBarry Smith .  man - manual page with additional information on option
994a264d7a6SBarry Smith .  list - the possible choices (one of these must be selected, anything else is invalid)
99553acd3b1SBarry Smith .  ntext - number of choices
9960fdccdaeSBarry Smith -  currentvalue - the current value; caller is responsible for setting this value correctly. Normally this is done with
9970fdccdaeSBarry Smith $                 PetscOptionsElist(..., obj->value,&value,&flg);
9980fdccdaeSBarry Smith $                 if (flg) {
9990fdccdaeSBarry Smith 
100053acd3b1SBarry Smith 
100153acd3b1SBarry Smith    Output Parameter:
100253acd3b1SBarry Smith +  value - the index of the value to return
100353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
100453acd3b1SBarry Smith 
100553acd3b1SBarry Smith    Level: intermediate
100653acd3b1SBarry Smith 
100753acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
100853acd3b1SBarry Smith 
1009a264d7a6SBarry Smith    See PetscOptionsFList() for when the choices are given in a PetscFunctionList()
101053acd3b1SBarry Smith 
101153acd3b1SBarry Smith    Concepts: options database^list
101253acd3b1SBarry Smith 
101353acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1014acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
101553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
101653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1017acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1018a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEnum()
101953acd3b1SBarry Smith @*/
1020*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsEList_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char ltext[],const char man[],const char *const *list,PetscInt ntext,const char currentvalue[],PetscInt *value,PetscBool  *set)
102153acd3b1SBarry Smith {
102253acd3b1SBarry Smith   PetscErrorCode ierr;
102353acd3b1SBarry Smith   PetscInt       i;
10241ae3d29cSBarry Smith   PetscOptions   amsopt;
102553acd3b1SBarry Smith 
102653acd3b1SBarry Smith   PetscFunctionBegin;
1027*1a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
1028*1a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,ltext,man,OPTION_ELIST,&amsopt);CHKERRQ(ierr);
102964facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
10300fdccdaeSBarry Smith     ierr = PetscStrdup(currentvalue ? currentvalue : "",(char**)&amsopt->data);CHKERRQ(ierr);
10311ae3d29cSBarry Smith     amsopt->list  = list;
10321ae3d29cSBarry Smith     amsopt->nlist = ntext;
10331ae3d29cSBarry Smith   }
1034*1a1499c8SBarry Smith   ierr = PetscOptionsGetEList(PetscOptionsObject->prefix,opt,list,ntext,value,set);CHKERRQ(ierr);
1035*1a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1036*1a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s> (choose one of)",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,currentvalue);CHKERRQ(ierr);
103753acd3b1SBarry Smith     for (i=0; i<ntext; i++) {
1038e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," %s",list[i]);CHKERRQ(ierr);
103953acd3b1SBarry Smith     }
1040e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm," (%s)\n",ManSection(man));CHKERRQ(ierr);
104153acd3b1SBarry Smith   }
104253acd3b1SBarry Smith   PetscFunctionReturn(0);
104353acd3b1SBarry Smith }
104453acd3b1SBarry Smith 
104553acd3b1SBarry Smith #undef __FUNCT__
1046e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupBegin_Private"
104753acd3b1SBarry Smith /*@C
1048acfcf0e5SJed Brown      PetscOptionsBoolGroupBegin - First in a series of logical queries on the options database for
1049d5649816SBarry Smith        which at most a single value can be true.
105053acd3b1SBarry Smith 
10513f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
105253acd3b1SBarry Smith 
105353acd3b1SBarry Smith    Input Parameters:
105453acd3b1SBarry Smith +  opt - option name
105553acd3b1SBarry Smith .  text - short string that describes the option
105653acd3b1SBarry Smith -  man - manual page with additional information on option
105753acd3b1SBarry Smith 
105853acd3b1SBarry Smith    Output Parameter:
105953acd3b1SBarry Smith .  flg - whether that option was set or not
106053acd3b1SBarry Smith 
106153acd3b1SBarry Smith    Level: intermediate
106253acd3b1SBarry Smith 
106353acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
106453acd3b1SBarry Smith 
1065acfcf0e5SJed Brown    Must be followed by 0 or more PetscOptionsBoolGroup()s and PetscOptionsBoolGroupEnd()
106653acd3b1SBarry Smith 
106753acd3b1SBarry Smith     Concepts: options database^logical group
106853acd3b1SBarry Smith 
106953acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1070acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
107153acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
107253acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1073acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1074a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
107553acd3b1SBarry Smith @*/
1076e55864a3SBarry Smith PetscErrorCode  PetscOptionsBoolGroupBegin_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
107753acd3b1SBarry Smith {
107853acd3b1SBarry Smith   PetscErrorCode ierr;
10791ae3d29cSBarry Smith   PetscOptions   amsopt;
108053acd3b1SBarry Smith 
108153acd3b1SBarry Smith   PetscFunctionBegin;
1082e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
108383355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1084ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1085a297a907SKarl Rupp 
1086ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
10871ae3d29cSBarry Smith   }
108868b16fdaSBarry Smith   *flg = PETSC_FALSE;
1089e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1090e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1091e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  Pick at most one of -------------\n");CHKERRQ(ierr);
1092e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
109353acd3b1SBarry Smith   }
109453acd3b1SBarry Smith   PetscFunctionReturn(0);
109553acd3b1SBarry Smith }
109653acd3b1SBarry Smith 
109753acd3b1SBarry Smith #undef __FUNCT__
1098e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroup_Private"
109953acd3b1SBarry Smith /*@C
1100acfcf0e5SJed Brown      PetscOptionsBoolGroup - One in a series of logical queries on the options database for
1101d5649816SBarry Smith        which at most a single value can be true.
110253acd3b1SBarry Smith 
11033f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
110453acd3b1SBarry Smith 
110553acd3b1SBarry Smith    Input Parameters:
110653acd3b1SBarry Smith +  opt - option name
110753acd3b1SBarry Smith .  text - short string that describes the option
110853acd3b1SBarry Smith -  man - manual page with additional information on option
110953acd3b1SBarry Smith 
111053acd3b1SBarry Smith    Output Parameter:
111153acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
111253acd3b1SBarry Smith 
111353acd3b1SBarry Smith    Level: intermediate
111453acd3b1SBarry Smith 
111553acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
111653acd3b1SBarry Smith 
1117acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin() and preceded a PetscOptionsBoolGroupEnd()
111853acd3b1SBarry Smith 
111953acd3b1SBarry Smith     Concepts: options database^logical group
112053acd3b1SBarry Smith 
112153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1122acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
112353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
112453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1125acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1126a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
112753acd3b1SBarry Smith @*/
1128e55864a3SBarry Smith PetscErrorCode  PetscOptionsBoolGroup_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
112953acd3b1SBarry Smith {
113053acd3b1SBarry Smith   PetscErrorCode ierr;
11311ae3d29cSBarry Smith   PetscOptions   amsopt;
113253acd3b1SBarry Smith 
113353acd3b1SBarry Smith   PetscFunctionBegin;
1134e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
113583355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1136ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1137a297a907SKarl Rupp 
1138ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11391ae3d29cSBarry Smith   }
114017326d04SJed Brown   *flg = PETSC_FALSE;
1141e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1142e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1143e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
114453acd3b1SBarry Smith   }
114553acd3b1SBarry Smith   PetscFunctionReturn(0);
114653acd3b1SBarry Smith }
114753acd3b1SBarry Smith 
114853acd3b1SBarry Smith #undef __FUNCT__
1149e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolGroupEnd_Private"
115053acd3b1SBarry Smith /*@C
1151acfcf0e5SJed Brown      PetscOptionsBoolGroupEnd - Last in a series of logical queries on the options database for
1152d5649816SBarry Smith        which at most a single value can be true.
115353acd3b1SBarry Smith 
11543f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
115553acd3b1SBarry Smith 
115653acd3b1SBarry Smith    Input Parameters:
115753acd3b1SBarry Smith +  opt - option name
115853acd3b1SBarry Smith .  text - short string that describes the option
115953acd3b1SBarry Smith -  man - manual page with additional information on option
116053acd3b1SBarry Smith 
116153acd3b1SBarry Smith    Output Parameter:
116253acd3b1SBarry Smith .  flg - PETSC_TRUE if found, else PETSC_FALSE
116353acd3b1SBarry Smith 
116453acd3b1SBarry Smith    Level: intermediate
116553acd3b1SBarry Smith 
116653acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
116753acd3b1SBarry Smith 
1168acfcf0e5SJed Brown    Must follow a PetscOptionsBoolGroupBegin()
116953acd3b1SBarry Smith 
117053acd3b1SBarry Smith     Concepts: options database^logical group
117153acd3b1SBarry Smith 
117253acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1173acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
117453acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
117553acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1176acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1177a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
117853acd3b1SBarry Smith @*/
1179e55864a3SBarry Smith PetscErrorCode  PetscOptionsBoolGroupEnd_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool  *flg)
118053acd3b1SBarry Smith {
118153acd3b1SBarry Smith   PetscErrorCode ierr;
11821ae3d29cSBarry Smith   PetscOptions   amsopt;
118353acd3b1SBarry Smith 
118453acd3b1SBarry Smith   PetscFunctionBegin;
1185e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
118683355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1187ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1188a297a907SKarl Rupp 
1189ace3abfcSBarry Smith     *(PetscBool*)amsopt->data = PETSC_FALSE;
11901ae3d29cSBarry Smith   }
119117326d04SJed Brown   *flg = PETSC_FALSE;
1192e55864a3SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,NULL);CHKERRQ(ierr);
1193e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1194e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"    -%s%s: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
119553acd3b1SBarry Smith   }
119653acd3b1SBarry Smith   PetscFunctionReturn(0);
119753acd3b1SBarry Smith }
119853acd3b1SBarry Smith 
119953acd3b1SBarry Smith #undef __FUNCT__
1200e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBool_Private"
120153acd3b1SBarry Smith /*@C
1202acfcf0e5SJed Brown    PetscOptionsBool - Determines if a particular option is in the database with a true or false
120353acd3b1SBarry Smith 
12043f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
120553acd3b1SBarry Smith 
120653acd3b1SBarry Smith    Input Parameters:
120753acd3b1SBarry Smith +  opt - option name
120853acd3b1SBarry Smith .  text - short string that describes the option
1209868c398cSBarry Smith .  man - manual page with additional information on option
121094ae4db5SBarry Smith -  currentvalue - the current value
121153acd3b1SBarry Smith 
121253acd3b1SBarry Smith    Output Parameter:
121353acd3b1SBarry Smith .  flg - PETSC_TRUE or PETSC_FALSE
121453acd3b1SBarry Smith .  set - PETSC_TRUE if found, else PETSC_FALSE
121553acd3b1SBarry Smith 
121653acd3b1SBarry Smith    Level: beginner
121753acd3b1SBarry Smith 
121853acd3b1SBarry Smith    Concepts: options database^logical
121953acd3b1SBarry Smith 
122053acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
122153acd3b1SBarry Smith 
122253acd3b1SBarry Smith .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
1223acfcf0e5SJed Brown           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1224acfcf0e5SJed Brown           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
122553acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
122653acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1227acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1228a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
122953acd3b1SBarry Smith @*/
1230*1a1499c8SBarry Smith PetscErrorCode  PetscOptionsBool_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool currentvalue,PetscBool  *flg,PetscBool  *set)
123153acd3b1SBarry Smith {
123253acd3b1SBarry Smith   PetscErrorCode ierr;
1233ace3abfcSBarry Smith   PetscBool      iset;
1234aee2cecaSBarry Smith   PetscOptions   amsopt;
123553acd3b1SBarry Smith 
123653acd3b1SBarry Smith   PetscFunctionBegin;
1237e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1238e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL,&amsopt);CHKERRQ(ierr);
1239ace3abfcSBarry Smith     ierr = PetscMalloc(sizeof(PetscBool),&amsopt->data);CHKERRQ(ierr);
1240a297a907SKarl Rupp 
124194ae4db5SBarry Smith     *(PetscBool*)amsopt->data = currentvalue;
1242af6d86caSBarry Smith   }
1243*1a1499c8SBarry Smith   ierr = PetscOptionsGetBool(PetscOptionsObject->prefix,opt,flg,&iset);CHKERRQ(ierr);
124453acd3b1SBarry Smith   if (set) *set = iset;
1245*1a1499c8SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
124694ae4db5SBarry Smith     const char *v = PetscBools[currentvalue];
1247*1a1499c8SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s: <%s> %s (%s)\n",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,v,text,ManSection(man));CHKERRQ(ierr);
124853acd3b1SBarry Smith   }
124953acd3b1SBarry Smith   PetscFunctionReturn(0);
125053acd3b1SBarry Smith }
125153acd3b1SBarry Smith 
125253acd3b1SBarry Smith #undef __FUNCT__
1253e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsRealArray_Private"
125453acd3b1SBarry Smith /*@C
125553acd3b1SBarry Smith    PetscOptionsRealArray - Gets an array of double values for a particular
125653acd3b1SBarry Smith    option in the database. The values must be separated with commas with
125753acd3b1SBarry Smith    no intervening spaces.
125853acd3b1SBarry Smith 
12593f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
126053acd3b1SBarry Smith 
126153acd3b1SBarry Smith    Input Parameters:
126253acd3b1SBarry Smith +  opt - the option one is seeking
126353acd3b1SBarry Smith .  text - short string describing option
126453acd3b1SBarry Smith .  man - manual page for option
126553acd3b1SBarry Smith -  nmax - maximum number of values
126653acd3b1SBarry Smith 
126753acd3b1SBarry Smith    Output Parameter:
126853acd3b1SBarry Smith +  value - location to copy values
126953acd3b1SBarry Smith .  nmax - actual number of values found
127053acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
127153acd3b1SBarry Smith 
127253acd3b1SBarry Smith    Level: beginner
127353acd3b1SBarry Smith 
127453acd3b1SBarry Smith    Notes:
127553acd3b1SBarry Smith    The user should pass in an array of doubles
127653acd3b1SBarry Smith 
127753acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
127853acd3b1SBarry Smith 
127953acd3b1SBarry Smith    Concepts: options database^array of strings
128053acd3b1SBarry Smith 
128153acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1282acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
128353acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
128453acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1285acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1286a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
128753acd3b1SBarry Smith @*/
1288e55864a3SBarry Smith PetscErrorCode  PetscOptionsRealArray_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscReal value[],PetscInt *n,PetscBool  *set)
128953acd3b1SBarry Smith {
129053acd3b1SBarry Smith   PetscErrorCode ierr;
129153acd3b1SBarry Smith   PetscInt       i;
1292e26ddf31SBarry Smith   PetscOptions   amsopt;
129353acd3b1SBarry Smith 
129453acd3b1SBarry Smith   PetscFunctionBegin;
1295e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1296e26ddf31SBarry Smith     PetscReal *vals;
1297e26ddf31SBarry Smith 
1298e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_REAL_ARRAY,&amsopt);CHKERRQ(ierr);
1299e55864a3SBarry Smith     ierr = PetscMalloc((*n)*sizeof(PetscReal),&amsopt->data);CHKERRQ(ierr);
1300e26ddf31SBarry Smith     vals = (PetscReal*)amsopt->data;
1301e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1302e26ddf31SBarry Smith     amsopt->arraylength = *n;
1303e26ddf31SBarry Smith   }
1304e55864a3SBarry Smith   ierr = PetscOptionsGetRealArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1305e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1306e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%G",PetscOptionsObject->prefix?PetscOptionsObject->prefix:"",opt+1,value[0]);CHKERRQ(ierr);
130753acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1308e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%g",(double)value[i]);CHKERRQ(ierr);
130953acd3b1SBarry Smith     }
1310e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
131153acd3b1SBarry Smith   }
131253acd3b1SBarry Smith   PetscFunctionReturn(0);
131353acd3b1SBarry Smith }
131453acd3b1SBarry Smith 
131553acd3b1SBarry Smith 
131653acd3b1SBarry Smith #undef __FUNCT__
1317e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsIntArray_Private"
131853acd3b1SBarry Smith /*@C
131953acd3b1SBarry Smith    PetscOptionsIntArray - Gets an array of integers for a particular
1320b32a342fSShri Abhyankar    option in the database.
132153acd3b1SBarry Smith 
13223f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
132353acd3b1SBarry Smith 
132453acd3b1SBarry Smith    Input Parameters:
132553acd3b1SBarry Smith +  opt - the option one is seeking
132653acd3b1SBarry Smith .  text - short string describing option
132753acd3b1SBarry Smith .  man - manual page for option
1328f8a50e2bSBarry Smith -  n - maximum number of values
132953acd3b1SBarry Smith 
133053acd3b1SBarry Smith    Output Parameter:
133153acd3b1SBarry Smith +  value - location to copy values
1332f8a50e2bSBarry Smith .  n - actual number of values found
133353acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
133453acd3b1SBarry Smith 
133553acd3b1SBarry Smith    Level: beginner
133653acd3b1SBarry Smith 
133753acd3b1SBarry Smith    Notes:
1338b32a342fSShri Abhyankar    The array can be passed as
1339b32a342fSShri Abhyankar    a comma seperated list:                                 0,1,2,3,4,5,6,7
13400fd488f5SShri Abhyankar    a range (start-end+1):                                  0-8
13410fd488f5SShri Abhyankar    a range with given increment (start-end+1:inc):         0-7:2
13420fd488f5SShri Abhyankar    a combination of values and ranges seperated by commas: 0,1-8,8-15:2
1343b32a342fSShri Abhyankar 
1344b32a342fSShri Abhyankar    There must be no intervening spaces between the values.
134553acd3b1SBarry Smith 
134653acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
134753acd3b1SBarry Smith 
1348b32a342fSShri Abhyankar    Concepts: options database^array of ints
134953acd3b1SBarry Smith 
135053acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1351acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
135253acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
135353acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1354acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1355a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsRealArray()
135653acd3b1SBarry Smith @*/
1357e55864a3SBarry Smith PetscErrorCode  PetscOptionsIntArray_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscInt value[],PetscInt *n,PetscBool  *set)
135853acd3b1SBarry Smith {
135953acd3b1SBarry Smith   PetscErrorCode ierr;
136053acd3b1SBarry Smith   PetscInt       i;
1361e26ddf31SBarry Smith   PetscOptions   amsopt;
136253acd3b1SBarry Smith 
136353acd3b1SBarry Smith   PetscFunctionBegin;
1364e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1365e26ddf31SBarry Smith     PetscInt *vals;
1366e26ddf31SBarry Smith 
1367e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_INT_ARRAY,&amsopt);CHKERRQ(ierr);
1368854ce69bSBarry Smith     ierr = PetscMalloc1(*n,(PetscInt**)&amsopt->data);CHKERRQ(ierr);
1369e26ddf31SBarry Smith     vals = (PetscInt*)amsopt->data;
1370e26ddf31SBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
1371e26ddf31SBarry Smith     amsopt->arraylength = *n;
1372e26ddf31SBarry Smith   }
1373e55864a3SBarry Smith   ierr = PetscOptionsGetIntArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1374e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1375e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
137653acd3b1SBarry Smith     for (i=1; i<*n; i++) {
1377e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
137853acd3b1SBarry Smith     }
1379e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
138053acd3b1SBarry Smith   }
138153acd3b1SBarry Smith   PetscFunctionReturn(0);
138253acd3b1SBarry Smith }
138353acd3b1SBarry Smith 
138453acd3b1SBarry Smith #undef __FUNCT__
1385e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsStringArray_Private"
138653acd3b1SBarry Smith /*@C
138753acd3b1SBarry Smith    PetscOptionsStringArray - Gets an array of string values for a particular
138853acd3b1SBarry Smith    option in the database. The values must be separated with commas with
138953acd3b1SBarry Smith    no intervening spaces.
139053acd3b1SBarry Smith 
13913f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
139253acd3b1SBarry Smith 
139353acd3b1SBarry Smith    Input Parameters:
139453acd3b1SBarry Smith +  opt - the option one is seeking
139553acd3b1SBarry Smith .  text - short string describing option
139653acd3b1SBarry Smith .  man - manual page for option
139753acd3b1SBarry Smith -  nmax - maximum number of strings
139853acd3b1SBarry Smith 
139953acd3b1SBarry Smith    Output Parameter:
140053acd3b1SBarry Smith +  value - location to copy strings
140153acd3b1SBarry Smith .  nmax - actual number of strings found
140253acd3b1SBarry Smith -  set - PETSC_TRUE if found, else PETSC_FALSE
140353acd3b1SBarry Smith 
140453acd3b1SBarry Smith    Level: beginner
140553acd3b1SBarry Smith 
140653acd3b1SBarry Smith    Notes:
140753acd3b1SBarry Smith    The user should pass in an array of pointers to char, to hold all the
140853acd3b1SBarry Smith    strings returned by this function.
140953acd3b1SBarry Smith 
141053acd3b1SBarry Smith    The user is responsible for deallocating the strings that are
141153acd3b1SBarry Smith    returned. The Fortran interface for this routine is not supported.
141253acd3b1SBarry Smith 
141353acd3b1SBarry Smith    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
141453acd3b1SBarry Smith 
141553acd3b1SBarry Smith    Concepts: options database^array of strings
141653acd3b1SBarry Smith 
141753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1418acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
141953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
142053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1421acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1422a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
142353acd3b1SBarry Smith @*/
1424e55864a3SBarry Smith PetscErrorCode  PetscOptionsStringArray_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],char *value[],PetscInt *nmax,PetscBool  *set)
142553acd3b1SBarry Smith {
142653acd3b1SBarry Smith   PetscErrorCode ierr;
14271ae3d29cSBarry Smith   PetscOptions   amsopt;
142853acd3b1SBarry Smith 
142953acd3b1SBarry Smith   PetscFunctionBegin;
1430e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1431e55864a3SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING_ARRAY,&amsopt);CHKERRQ(ierr);
1432854ce69bSBarry Smith     ierr = PetscMalloc1(*nmax,(char**)&amsopt->data);CHKERRQ(ierr);
1433a297a907SKarl Rupp 
14341ae3d29cSBarry Smith     amsopt->arraylength = *nmax;
14351ae3d29cSBarry Smith   }
1436e55864a3SBarry Smith   ierr = PetscOptionsGetStringArray(PetscOptionsObject->prefix,opt,value,nmax,set);CHKERRQ(ierr);
1437e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1438e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <string1,string2,...>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,text,ManSection(man));CHKERRQ(ierr);
143953acd3b1SBarry Smith   }
144053acd3b1SBarry Smith   PetscFunctionReturn(0);
144153acd3b1SBarry Smith }
144253acd3b1SBarry Smith 
1443e2446a98SMatthew Knepley #undef __FUNCT__
1444e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsBoolArray_Private"
1445e2446a98SMatthew Knepley /*@C
1446acfcf0e5SJed Brown    PetscOptionsBoolArray - Gets an array of logical values (true or false) for a particular
1447e2446a98SMatthew Knepley    option in the database. The values must be separated with commas with
1448e2446a98SMatthew Knepley    no intervening spaces.
1449e2446a98SMatthew Knepley 
14503f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
1451e2446a98SMatthew Knepley 
1452e2446a98SMatthew Knepley    Input Parameters:
1453e2446a98SMatthew Knepley +  opt - the option one is seeking
1454e2446a98SMatthew Knepley .  text - short string describing option
1455e2446a98SMatthew Knepley .  man - manual page for option
1456e2446a98SMatthew Knepley -  nmax - maximum number of values
1457e2446a98SMatthew Knepley 
1458e2446a98SMatthew Knepley    Output Parameter:
1459e2446a98SMatthew Knepley +  value - location to copy values
1460e2446a98SMatthew Knepley .  nmax - actual number of values found
1461e2446a98SMatthew Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
1462e2446a98SMatthew Knepley 
1463e2446a98SMatthew Knepley    Level: beginner
1464e2446a98SMatthew Knepley 
1465e2446a98SMatthew Knepley    Notes:
1466e2446a98SMatthew Knepley    The user should pass in an array of doubles
1467e2446a98SMatthew Knepley 
1468e2446a98SMatthew Knepley    Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
1469e2446a98SMatthew Knepley 
1470e2446a98SMatthew Knepley    Concepts: options database^array of strings
1471e2446a98SMatthew Knepley 
1472e2446a98SMatthew Knepley .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1473acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
1474e2446a98SMatthew Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1475e2446a98SMatthew Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1476acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1477a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
1478e2446a98SMatthew Knepley @*/
1479e55864a3SBarry Smith PetscErrorCode  PetscOptionsBoolArray_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscBool value[],PetscInt *n,PetscBool *set)
1480e2446a98SMatthew Knepley {
1481e2446a98SMatthew Knepley   PetscErrorCode ierr;
1482e2446a98SMatthew Knepley   PetscInt       i;
14831ae3d29cSBarry Smith   PetscOptions   amsopt;
1484e2446a98SMatthew Knepley 
1485e2446a98SMatthew Knepley   PetscFunctionBegin;
1486e55864a3SBarry Smith   if (!PetscOptionsObject->count) {
1487ace3abfcSBarry Smith     PetscBool *vals;
14881ae3d29cSBarry Smith 
148983355fc5SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_BOOL_ARRAY,&amsopt);CHKERRQ(ierr);
1490*1a1499c8SBarry Smith     ierr = PetscMalloc1(*n,(PetscBool**)&amsopt->data);CHKERRQ(ierr);
1491ace3abfcSBarry Smith     vals = (PetscBool*)amsopt->data;
14921ae3d29cSBarry Smith     for (i=0; i<*n; i++) vals[i] = value[i];
14931ae3d29cSBarry Smith     amsopt->arraylength = *n;
14941ae3d29cSBarry Smith   }
1495e55864a3SBarry Smith   ierr = PetscOptionsGetBoolArray(PetscOptionsObject->prefix,opt,value,n,set);CHKERRQ(ierr);
1496e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1497e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%d",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,value[0]);CHKERRQ(ierr);
1498e2446a98SMatthew Knepley     for (i=1; i<*n; i++) {
1499e55864a3SBarry Smith       ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,",%d",value[i]);CHKERRQ(ierr);
1500e2446a98SMatthew Knepley     }
1501e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,">: %s (%s)\n",text,ManSection(man));CHKERRQ(ierr);
1502e2446a98SMatthew Knepley   }
1503e2446a98SMatthew Knepley   PetscFunctionReturn(0);
1504e2446a98SMatthew Knepley }
1505e2446a98SMatthew Knepley 
15068cc676e6SMatthew G Knepley #undef __FUNCT__
1507e55864a3SBarry Smith #define __FUNCT__ "PetscOptionsViewer_Private"
15088cc676e6SMatthew G Knepley /*@C
1509d1da0b69SBarry Smith    PetscOptionsViewer - Gets a viewer appropriate for the type indicated by the user
15108cc676e6SMatthew G Knepley 
15118cc676e6SMatthew G Knepley    Logically Collective on the communicator passed in PetscOptionsBegin()
15128cc676e6SMatthew G Knepley 
15138cc676e6SMatthew G Knepley    Input Parameters:
15148cc676e6SMatthew G Knepley +  opt - option name
15158cc676e6SMatthew G Knepley .  text - short string that describes the option
15168cc676e6SMatthew G Knepley -  man - manual page with additional information on option
15178cc676e6SMatthew G Knepley 
15188cc676e6SMatthew G Knepley    Output Parameter:
15198cc676e6SMatthew G Knepley +  viewer - the viewer
15208cc676e6SMatthew G Knepley -  set - PETSC_TRUE if found, else PETSC_FALSE
15218cc676e6SMatthew G Knepley 
15228cc676e6SMatthew G Knepley    Level: beginner
15238cc676e6SMatthew G Knepley 
15248cc676e6SMatthew G Knepley    Concepts: options database^has int
15258cc676e6SMatthew G Knepley 
15268cc676e6SMatthew G Knepley    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
15278cc676e6SMatthew G Knepley 
15285a7113b9SPatrick Sanan    See PetscOptionsGetViewer() for the format of the supplied viewer and its options
15298cc676e6SMatthew G Knepley 
15308cc676e6SMatthew G Knepley .seealso: PetscOptionsGetViewer(), PetscOptionsHasName(), PetscOptionsGetString(), PetscOptionsGetInt(),
15318cc676e6SMatthew G Knepley           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
15328cc676e6SMatthew G Knepley           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
15338cc676e6SMatthew G Knepley           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
15348cc676e6SMatthew G Knepley           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
15358cc676e6SMatthew G Knepley           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1536a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
15378cc676e6SMatthew G Knepley @*/
1538e55864a3SBarry Smith PetscErrorCode  PetscOptionsViewer_Private(PetscOptionsObjectType *PetscOptionsObject,const char opt[],const char text[],const char man[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
15398cc676e6SMatthew G Knepley {
15408cc676e6SMatthew G Knepley   PetscErrorCode ierr;
15418cc676e6SMatthew G Knepley   PetscOptions   amsopt;
15428cc676e6SMatthew G Knepley 
15438cc676e6SMatthew G Knepley   PetscFunctionBegin;
1544*1a1499c8SBarry Smith   if (!PetscOptionsObject->count) {
1545*1a1499c8SBarry Smith     ierr = PetscOptionsCreate_Private(PetscOptionsObject,opt,text,man,OPTION_STRING,&amsopt);CHKERRQ(ierr);
154664facd6cSBarry Smith     /* must use system malloc since SAWs may free this */
15475b02f95dSBarry Smith     ierr = PetscStrdup("",(char**)&amsopt->data);CHKERRQ(ierr);
15488cc676e6SMatthew G Knepley   }
1549e55864a3SBarry Smith   ierr = PetscOptionsGetViewer(PetscOptionsObject->comm,PetscOptionsObject->prefix,opt,viewer,format,set);CHKERRQ(ierr);
1550e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1551e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  -%s%s <%s>: %s (%s)\n",PetscOptionsObject->prefix ? PetscOptionsObject->prefix : "",opt+1,"",text,ManSection(man));CHKERRQ(ierr);
15528cc676e6SMatthew G Knepley   }
15538cc676e6SMatthew G Knepley   PetscFunctionReturn(0);
15548cc676e6SMatthew G Knepley }
15558cc676e6SMatthew G Knepley 
155653acd3b1SBarry Smith 
155753acd3b1SBarry Smith #undef __FUNCT__
155853acd3b1SBarry Smith #define __FUNCT__ "PetscOptionsHead"
155953acd3b1SBarry Smith /*@C
1560b52f573bSBarry Smith      PetscOptionsHead - Puts a heading before listing any more published options. Used, for example,
156153acd3b1SBarry Smith             in KSPSetFromOptions_GMRES().
156253acd3b1SBarry Smith 
15633f9fe445SBarry Smith    Logically Collective on the communicator passed in PetscOptionsBegin()
156453acd3b1SBarry Smith 
156553acd3b1SBarry Smith    Input Parameter:
156653acd3b1SBarry Smith .   head - the heading text
156753acd3b1SBarry Smith 
156853acd3b1SBarry Smith 
156953acd3b1SBarry Smith    Level: intermediate
157053acd3b1SBarry Smith 
157153acd3b1SBarry Smith    Notes: Must be between a PetscOptionsBegin() and a PetscOptionsEnd()
157253acd3b1SBarry Smith 
1573b52f573bSBarry Smith           Can be followed by a call to PetscOptionsTail() in the same function.
157453acd3b1SBarry Smith 
157553acd3b1SBarry Smith    Concepts: options database^subheading
157653acd3b1SBarry Smith 
157753acd3b1SBarry Smith .seealso: PetscOptionsGetInt(), PetscOptionsGetReal(),
1578acfcf0e5SJed Brown            PetscOptionsHasName(), PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool(),
157953acd3b1SBarry Smith           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
158053acd3b1SBarry Smith           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1581acfcf0e5SJed Brown           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1582a264d7a6SBarry Smith           PetscOptionsFList(), PetscOptionsEList()
158353acd3b1SBarry Smith @*/
1584e55864a3SBarry Smith PetscErrorCode  PetscOptionsHead(PetscOptionsObjectType *PetscOptionsObject,const char head[])
158553acd3b1SBarry Smith {
158653acd3b1SBarry Smith   PetscErrorCode ierr;
158753acd3b1SBarry Smith 
158853acd3b1SBarry Smith   PetscFunctionBegin;
1589e55864a3SBarry Smith   if (PetscOptionsObject->printhelp && PetscOptionsObject->count == 1 && !PetscOptionsObject->alreadyprinted) {
1590e55864a3SBarry Smith     ierr = (*PetscHelpPrintf)(PetscOptionsObject->comm,"  %s\n",head);CHKERRQ(ierr);
159153acd3b1SBarry Smith   }
159253acd3b1SBarry Smith   PetscFunctionReturn(0);
159353acd3b1SBarry Smith }
159453acd3b1SBarry Smith 
159553acd3b1SBarry Smith 
159653acd3b1SBarry Smith 
159753acd3b1SBarry Smith 
159853acd3b1SBarry Smith 
159953acd3b1SBarry Smith 
1600